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"),
224 fUseTrackingErrorsForAlignment(""),
225 fGAliceFileName(gAliceFilename),
230 fNumberOfEventsPerFile(1),
232 fLoadAlignFromCDB(kTRUE),
233 fLoadAlignData("ALL"),
241 fParentRawReader(NULL),
246 fDiamondProfile(NULL),
247 fDiamondProfileTPC(NULL),
248 fMeanVertexConstraint(kTRUE),
252 fAlignObjArray(NULL),
255 fInitCDBCalled(kFALSE),
256 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;
287 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
291 //_____________________________________________________________________________
292 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
294 fUniformField(rec.fUniformField),
295 fForcedFieldMap(NULL),
296 fRunVertexFinder(rec.fRunVertexFinder),
297 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
298 fRunHLTTracking(rec.fRunHLTTracking),
299 fRunMuonTracking(rec.fRunMuonTracking),
300 fRunV0Finder(rec.fRunV0Finder),
301 fRunCascadeFinder(rec.fRunCascadeFinder),
302 fStopOnError(rec.fStopOnError),
303 fWriteAlignmentData(rec.fWriteAlignmentData),
304 fWriteESDfriend(rec.fWriteESDfriend),
305 fFillTriggerESD(rec.fFillTriggerESD),
307 fCleanESD(rec.fCleanESD),
308 fV0DCAmax(rec.fV0DCAmax),
309 fV0CsPmin(rec.fV0CsPmin),
313 fRunLocalReconstruction(rec.fRunLocalReconstruction),
314 fRunTracking(rec.fRunTracking),
315 fFillESD(rec.fFillESD),
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 fInLoopQA(rec.fInLoopQA),
355 fSameQACycle(rec.fSameQACycle),
356 fRunPlaneEff(rec.fRunPlaneEff),
365 fIsNewRunLoader(rec.fIsNewRunLoader),
371 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
372 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
374 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
375 fReconstructor[iDet] = NULL;
376 fLoader[iDet] = NULL;
377 fTracker[iDet] = NULL;
378 fQACycles[iDet] = rec.fQACycles[iDet];
380 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
381 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
386 //_____________________________________________________________________________
387 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
389 // assignment operator
390 // Used in PROOF mode
391 // Be very careful while modifing it!
392 // Simple rules to follow:
393 // for persistent data members - use their assignment operators
394 // for non-persistent ones - do nothing or take the default values from constructor
395 // TSelector members should not be touched
396 if(&rec == this) return *this;
398 fUniformField = rec.fUniformField;
399 fForcedFieldMap = NULL;
400 fRunVertexFinder = rec.fRunVertexFinder;
401 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
402 fRunHLTTracking = rec.fRunHLTTracking;
403 fRunMuonTracking = rec.fRunMuonTracking;
404 fRunV0Finder = rec.fRunV0Finder;
405 fRunCascadeFinder = rec.fRunCascadeFinder;
406 fStopOnError = rec.fStopOnError;
407 fWriteAlignmentData = rec.fWriteAlignmentData;
408 fWriteESDfriend = rec.fWriteESDfriend;
409 fFillTriggerESD = rec.fFillTriggerESD;
411 fCleanESD = rec.fCleanESD;
412 fV0DCAmax = rec.fV0DCAmax;
413 fV0CsPmin = rec.fV0CsPmin;
417 fRunLocalReconstruction = rec.fRunLocalReconstruction;
418 fRunTracking = rec.fRunTracking;
419 fFillESD = rec.fFillESD;
420 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
421 fGAliceFileName = rec.fGAliceFileName;
422 fRawInput = rec.fRawInput;
423 fEquipIdMap = rec.fEquipIdMap;
424 fFirstEvent = rec.fFirstEvent;
425 fLastEvent = rec.fLastEvent;
426 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
428 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
429 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
432 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
433 fLoadAlignData = rec.fLoadAlignData;
434 fESDPar = rec.fESDPar;
435 fUseHLTData = rec.fUseHLTData;
437 delete fRunInfo; fRunInfo = NULL;
438 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
440 fEventInfo = rec.fEventInfo;
444 fParentRawReader = NULL;
446 fRecoParam = rec.fRecoParam;
448 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
449 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
450 delete fLoader[iDet]; fLoader[iDet] = NULL;
451 delete fTracker[iDet]; fTracker[iDet] = NULL;
452 fQACycles[iDet] = rec.fQACycles[iDet];
456 delete fDiamondProfile; fDiamondProfile = NULL;
457 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
458 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
459 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
460 fMeanVertexConstraint = rec.fMeanVertexConstraint;
462 delete fGRPData; fGRPData = NULL;
463 if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
465 delete fAlignObjArray; fAlignObjArray = NULL;
468 fSpecCDBUri.Delete();
469 fInitCDBCalled = rec.fInitCDBCalled;
470 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
471 fQADetectors = rec.fQADetectors;
473 fQATasks = rec.fQATasks;
475 fRunGlobalQA = rec.fRunGlobalQA;
476 fInLoopQA = rec.fInLoopQA;
477 fSameQACycle = rec.fSameQACycle;
478 fRunPlaneEff = rec.fRunPlaneEff;
487 fIsNewRunLoader = rec.fIsNewRunLoader;
494 //_____________________________________________________________________________
495 AliReconstruction::~AliReconstruction()
500 delete fForcedFieldMap;
502 if (fAlignObjArray) {
503 fAlignObjArray->Delete();
504 delete fAlignObjArray;
506 fSpecCDBUri.Delete();
508 AliCodeTimer::Instance()->Print();
511 //_____________________________________________________________________________
512 void AliReconstruction::InitCDB()
514 // activate a default CDB storage
515 // First check if we have any CDB storage set, because it is used
516 // to retrieve the calibration and alignment constants
517 AliCodeTimerAuto("");
519 if (fInitCDBCalled) return;
520 fInitCDBCalled = kTRUE;
522 AliCDBManager* man = AliCDBManager::Instance();
523 if (man->IsDefaultStorageSet())
525 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
526 AliWarning("Default CDB storage has been already set !");
527 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
528 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
529 fCDBUri = man->GetDefaultStorage()->GetURI();
532 if (fCDBUri.Length() > 0)
534 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
535 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
536 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
538 fCDBUri="local://$ALICE_ROOT";
539 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
540 AliWarning("Default CDB storage not yet set !!!!");
541 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
542 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
545 man->SetDefaultStorage(fCDBUri);
548 // Now activate the detector specific CDB storage locations
549 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
550 TObject* obj = fSpecCDBUri[i];
552 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
553 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
554 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
555 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
557 AliSysInfo::AddStamp("InitCDB");
560 //_____________________________________________________________________________
561 void AliReconstruction::SetDefaultStorage(const char* uri) {
562 // Store the desired default CDB storage location
563 // Activate it later within the Run() method
569 //_____________________________________________________________________________
570 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
571 // Store a detector-specific CDB storage location
572 // Activate it later within the Run() method
574 AliCDBPath aPath(calibType);
575 if(!aPath.IsValid()){
576 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
577 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
578 if(!strcmp(calibType, fgkDetectorName[iDet])) {
579 aPath.SetPath(Form("%s/*", calibType));
580 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
584 if(!aPath.IsValid()){
585 AliError(Form("Not a valid path or detector: %s", calibType));
590 // // check that calibType refers to a "valid" detector name
591 // Bool_t isDetector = kFALSE;
592 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
593 // TString detName = fgkDetectorName[iDet];
594 // if(aPath.GetLevel0() == detName) {
595 // isDetector = kTRUE;
601 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
605 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
606 if (obj) fSpecCDBUri.Remove(obj);
607 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
611 //_____________________________________________________________________________
612 Bool_t AliReconstruction::SetRunNumberFromData()
614 // The method is called in Run() in order
615 // to set a correct run number.
616 // In case of raw data reconstruction the
617 // run number is taken from the raw data header
619 if (fSetRunNumberFromDataCalled) return kTRUE;
620 fSetRunNumberFromDataCalled = kTRUE;
622 AliCDBManager* man = AliCDBManager::Instance();
625 if(fRawReader->NextEvent()) {
626 if(man->GetRun() > 0) {
627 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
629 man->SetRun(fRawReader->GetRunNumber());
630 fRawReader->RewindEvents();
633 if(man->GetRun() > 0) {
634 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
637 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
643 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
645 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
650 // read run number from gAlice
651 if(rl->GetHeader()) {
652 man->SetRun(rl->GetHeader()->GetRun());
657 AliError("Neither run-loader header nor RawReader objects are found !");
669 //_____________________________________________________________________________
670 void AliReconstruction::SetCDBLock() {
671 // Set CDB lock: from now on it is forbidden to reset the run number
672 // or the default storage or to activate any further storage!
674 AliCDBManager::Instance()->SetLock(1);
677 //_____________________________________________________________________________
678 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
680 // Read the alignment objects from CDB.
681 // Each detector is supposed to have the
682 // alignment objects in DET/Align/Data CDB path.
683 // All the detector objects are then collected,
684 // sorted by geometry level (starting from ALIC) and
685 // then applied to the TGeo geometry.
686 // Finally an overlaps check is performed.
688 // Load alignment data from CDB and fill fAlignObjArray
689 if(fLoadAlignFromCDB){
691 TString detStr = detectors;
692 TString loadAlObjsListOfDets = "";
694 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
695 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
696 loadAlObjsListOfDets += fgkDetectorName[iDet];
697 loadAlObjsListOfDets += " ";
698 } // end loop over detectors
699 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
700 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
702 // Check if the array with alignment objects was
703 // provided by the user. If yes, apply the objects
704 // to the present TGeo geometry
705 if (fAlignObjArray) {
706 if (gGeoManager && gGeoManager->IsClosed()) {
707 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
708 AliError("The misalignment of one or more volumes failed!"
709 "Compare the list of simulated detectors and the list of detector alignment data!");
714 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
720 if (fAlignObjArray) {
721 fAlignObjArray->Delete();
722 delete fAlignObjArray; fAlignObjArray=NULL;
728 //_____________________________________________________________________________
729 void AliReconstruction::SetGAliceFile(const char* fileName)
731 // set the name of the galice file
733 fGAliceFileName = fileName;
736 //_____________________________________________________________________________
737 void AliReconstruction::SetInput(const char* input)
739 // In case the input string starts with 'mem://', we run in an online mode
740 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
741 // file is assumed. One can give as an input:
742 // mem://: - events taken from DAQ monitoring libs online
744 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
745 if (input) fRawInput = input;
748 //_____________________________________________________________________________
749 void AliReconstruction::SetOption(const char* detector, const char* option)
751 // set options for the reconstruction of a detector
753 TObject* obj = fOptions.FindObject(detector);
754 if (obj) fOptions.Remove(obj);
755 fOptions.Add(new TNamed(detector, option));
758 //_____________________________________________________________________________
759 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
761 // Set custom reconstruction parameters for a given detector
762 // Single set of parameters for all the events
763 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
764 if(!strcmp(detector, fgkDetectorName[iDet])) {
766 fRecoParam.AddDetRecoParam(iDet,par);
773 //_____________________________________________________________________________
774 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
775 //------------------------------------------------
776 // The magnetic field map, defined externally...
777 // L3 current 30000 A -> 0.5 T
778 // L3 current 12000 A -> 0.2 T
779 // dipole current 6000 A
780 // The polarities must be the same
781 //------------------------------------------------
782 const Float_t l3NominalCurrent1=30000.; // (A)
783 const Float_t l3NominalCurrent2=12000.; // (A)
784 const Float_t diNominalCurrent =6000. ; // (A)
786 const Float_t tolerance=0.03; // relative current tolerance
787 const Float_t zero=77.; // "zero" current (A)
790 Bool_t dipoleON=kFALSE;
792 TString s=(factor < 0) ? "L3: -" : "L3: +";
794 l3Current = TMath::Abs(l3Current);
795 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
796 map=AliMagWrapCheb::k5kG;
799 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
800 map=AliMagWrapCheb::k2kG;
803 if (l3Current < zero) {
804 map=AliMagWrapCheb::k2kG;
806 factor=0.; // in fact, this is a global factor...
807 fUniformField=kTRUE; // track with the uniform (zero) B field
809 AliError(Form("Wrong L3 current (%f A)!",l3Current));
813 diCurrent = TMath::Abs(diCurrent);
814 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
815 // 3% current tolerance...
819 if (diCurrent < zero) { // some small current..
823 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
827 delete fForcedFieldMap;
829 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
831 fForcedFieldMap->Print();
833 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
839 Bool_t AliReconstruction::InitGRP() {
840 //------------------------------------
841 // Initialization of the GRP entry
842 //------------------------------------
843 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
846 fGRPData = dynamic_cast<TMap*>(entry->GetObject());
848 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
852 AliError("No GRP entry found in OCDB!");
856 TObjString *lhcState=
857 dynamic_cast<TObjString*>(fGRPData->GetValue("fLHCState"));
859 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
862 TObjString *beamType=
863 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamType"));
865 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
868 TObjString *beamEnergyStr=
869 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamEnergy"));
870 if (!beamEnergyStr) {
871 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
875 dynamic_cast<TObjString*>(fGRPData->GetValue("fRunType"));
877 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
880 TObjString *activeDetectors=
881 dynamic_cast<TObjString*>(fGRPData->GetValue("fDetectorMask"));
882 if (!activeDetectors) {
883 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
886 fRunInfo = new AliRunInfo(lhcState ? lhcState->GetString().Data() : "UNKNOWN",
887 beamType ? beamType->GetString().Data() : "UNKNOWN",
888 beamEnergyStr ? beamEnergyStr->GetString().Atof() : 0,
889 runType ? runType->GetString().Data() : "UNKNOWN",
890 activeDetectors ? activeDetectors->GetString().Atoi() : 1074790399);
892 // Process the list of active detectors
893 if (activeDetectors && activeDetectors->GetString().IsDigit()) {
894 UInt_t detMask = activeDetectors->GetString().Atoi();
895 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
896 fRunTracking = MatchDetectorList(fRunTracking,detMask);
897 fFillESD = MatchDetectorList(fFillESD,detMask);
898 fQADetectors = MatchDetectorList(fQADetectors,detMask);
901 AliInfo("===================================================================================");
902 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
903 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
904 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
905 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
906 AliInfo("===================================================================================");
908 //*** Dealing with the magnetic field map
909 if (AliTracker::GetFieldMap()) {
910 AliInfo("Running with the externally set B field !");
912 // Construct the field map out of the information retrieved from GRP.
917 TObjString *l3Current=
918 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
920 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
923 TObjString *l3Polarity=
924 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
926 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
931 TObjString *diCurrent=
932 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
934 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
937 TObjString *diPolarity=
938 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
940 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
945 Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
946 Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
947 Float_t l3Pol=atof(l3Polarity->GetName());
949 if (l3Pol != 0.) factor=-1.;
952 if (!SetFieldMap(l3Cur, diCur, factor)) {
953 AliFatal("Failed to creat a B field map ! Exiting...");
955 AliInfo("Running with the B field constructed out of GRP !");
958 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
963 //*** Get the diamond profile from OCDB
964 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
966 if (fMeanVertexConstraint)
967 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
969 AliError("No diamond profile found in OCDB!");
972 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
974 if (fMeanVertexConstraint)
975 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
977 AliError("No diamond profile found in OCDB!");
983 //_____________________________________________________________________________
984 Bool_t AliReconstruction::LoadCDB()
986 AliCodeTimerAuto("");
988 AliCDBManager::Instance()->Get("GRP/CTP/Config");
990 TString detStr = fRunLocalReconstruction;
991 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
992 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
993 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
998 //_____________________________________________________________________________
999 Bool_t AliReconstruction::Run(const char* input)
1002 AliCodeTimerAuto("");
1005 if (GetAbort() != TSelector::kContinue) return kFALSE;
1007 TChain *chain = NULL;
1008 if (fRawReader && (chain = fRawReader->GetChain())) {
1011 gProof->AddInput(this);
1013 outputFile.SetProtocol("root",kTRUE);
1014 outputFile.SetHost(gSystem->HostName());
1015 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1016 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1017 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1019 chain->Process("AliReconstruction");
1022 chain->Process(this);
1027 if (GetAbort() != TSelector::kContinue) return kFALSE;
1029 if (GetAbort() != TSelector::kContinue) return kFALSE;
1030 //******* The loop over events
1032 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1033 (fRawReader && fRawReader->NextEvent())) {
1034 if (!ProcessEvent(iEvent)) {
1035 Abort("ProcessEvent",TSelector::kAbortFile);
1041 if (GetAbort() != TSelector::kContinue) return kFALSE;
1043 if (GetAbort() != TSelector::kContinue) return kFALSE;
1049 //_____________________________________________________________________________
1050 void AliReconstruction::InitRawReader(const char* input)
1052 AliCodeTimerAuto("");
1054 // Init raw-reader and
1055 // set the input in case of raw data
1056 if (input) fRawInput = input;
1057 fRawReader = AliRawReader::Create(fRawInput.Data());
1059 AliInfo("Reconstruction will run over digits");
1061 if (!fEquipIdMap.IsNull() && fRawReader)
1062 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1064 if (!fUseHLTData.IsNull()) {
1065 // create the RawReaderHLT which performs redirection of HLT input data for
1066 // the specified detectors
1067 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1069 fParentRawReader=fRawReader;
1070 fRawReader=pRawReader;
1072 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1075 AliSysInfo::AddStamp("CreateRawReader");
1078 //_____________________________________________________________________________
1079 void AliReconstruction::InitRun(const char* input)
1081 // Initialization of raw-reader,
1082 // run number, CDB etc.
1083 AliCodeTimerAuto("");
1084 AliSysInfo::AddStamp("Start");
1086 // Initialize raw-reader if any
1087 InitRawReader(input);
1089 // Initialize the CDB storage
1092 // Set run number in CDBManager (if it is not already set by the user)
1093 if (!SetRunNumberFromData()) {
1094 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1098 // Set CDB lock: from now on it is forbidden to reset the run number
1099 // or the default storage or to activate any further storage!
1104 //_____________________________________________________________________________
1105 void AliReconstruction::Begin(TTree *)
1107 // Initialize AlReconstruction before
1108 // going into the event loop
1109 // Should follow the TSelector convention
1110 // i.e. initialize only the object on the client side
1111 AliCodeTimerAuto("");
1113 AliReconstruction *reco = NULL;
1115 if (reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1118 AliSysInfo::AddStamp("ReadInputInBegin");
1121 // Import ideal TGeo geometry and apply misalignment
1123 TString geom(gSystem->DirName(fGAliceFileName));
1124 geom += "/geometry.root";
1125 AliGeomManager::LoadGeometry(geom.Data());
1127 Abort("LoadGeometry", TSelector::kAbortProcess);
1130 AliSysInfo::AddStamp("LoadGeom");
1131 TString detsToCheck=fRunLocalReconstruction;
1132 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1133 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1136 AliSysInfo::AddStamp("CheckGeom");
1139 if (!MisalignGeometry(fLoadAlignData)) {
1140 Abort("MisalignGeometry", TSelector::kAbortProcess);
1143 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1144 AliSysInfo::AddStamp("MisalignGeom");
1147 Abort("InitGRP", TSelector::kAbortProcess);
1150 AliSysInfo::AddStamp("InitGRP");
1153 Abort("LoadCDB", TSelector::kAbortProcess);
1156 AliSysInfo::AddStamp("LoadCDB");
1158 // Read the reconstruction parameters from OCDB
1159 if (!InitRecoParams()) {
1160 AliWarning("Not all detectors have correct RecoParam objects initialized");
1162 AliSysInfo::AddStamp("InitRecoParams");
1165 if (reco) *reco = *this;
1166 fInput->Add(gGeoManager);
1168 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1169 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1170 fInput->Add((AliMagF*)AliTracker::GetFieldMap());
1175 //_____________________________________________________________________________
1176 void AliReconstruction::SlaveBegin(TTree*)
1178 // Initialization related to run-loader,
1179 // vertexer, trackers, recontructors
1180 // In proof mode it is executed on the slave
1181 AliCodeTimerAuto("");
1183 TProofOutputFile *outProofFile = NULL;
1185 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1188 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1190 AliGeomManager::SetGeometry(tgeo);
1192 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1193 Int_t runNumber = -1;
1194 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1195 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1196 man->SetCacheFlag(kTRUE);
1197 man->SetLock(kTRUE);
1201 if (AliMagF *map = (AliMagF*)fInput->FindObject("Maps")) {
1202 AliTracker::SetFieldMap(map,fUniformField);
1204 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1205 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1206 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1207 fOutput->Add(outProofFile);
1209 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1212 // get the run loader
1213 if (!InitRunLoader()) {
1214 Abort("InitRunLoader", TSelector::kAbortProcess);
1217 AliSysInfo::AddStamp("LoadLoader");
1219 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1220 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
1223 if (fRunVertexFinder && !CreateVertexer()) {
1224 Abort("CreateVertexer", TSelector::kAbortProcess);
1227 AliSysInfo::AddStamp("CreateVertexer");
1230 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1231 Abort("CreateTrackers", TSelector::kAbortProcess);
1234 AliSysInfo::AddStamp("CreateTrackers");
1236 // create the ESD output file and tree
1237 if (!outProofFile) {
1238 ffile = TFile::Open("AliESDs.root", "RECREATE");
1239 ffile->SetCompressionLevel(2);
1240 if (!ffile->IsOpen()) {
1241 Abort("OpenESDFile", TSelector::kAbortProcess);
1246 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1247 Abort(Form("Problems opening output PROOF file: %s/%s",
1248 outProofFile->GetDir(), outProofFile->GetFileName()),
1249 TSelector::kAbortProcess);
1254 ftree = new TTree("esdTree", "Tree with ESD objects");
1255 fesd = new AliESDEvent();
1256 fesd->CreateStdContent();
1257 fesd->WriteToTree(ftree);
1259 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1260 fhltesd = new AliESDEvent();
1261 fhltesd->CreateStdContent();
1262 fhltesd->WriteToTree(fhlttree);
1265 if (fWriteESDfriend) {
1266 fesdf = new AliESDfriend();
1267 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1268 br->SetFile("AliESDfriends.root");
1269 fesd->AddObject(fesdf);
1272 ProcInfo_t ProcInfo;
1273 gSystem->GetProcInfo(&ProcInfo);
1274 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1277 fQASteer = new AliQADataMakerSteer("rec") ;
1278 fQASteer->SetActiveDetectors(fQADetectors) ;
1279 fQASteer->SetTasks(fQATasks) ;
1282 if (fRunQA && fRawReader && fQATasks.Contains(Form("%d", AliQA::kRAWS))) {
1283 fQASteer->Run(fQADetectors, fRawReader) ;
1284 fSameQACycle = kTRUE ;
1288 //Initialize the QA and start of cycle for out-of-loop QA
1290 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam, fSameQACycle, !fInLoopQA) ;
1294 fSameQACycle = kFALSE;
1295 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1296 AliInfo(Form("Initializing the global QA data maker"));
1297 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1298 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
1299 AliTracker::SetResidualsArray(arr);
1301 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1302 fSameQACycle = kTRUE;
1305 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1306 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
1308 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1309 fSameQACycle = kTRUE;
1314 //Initialize the Plane Efficiency framework
1315 if (fRunPlaneEff && !InitPlaneEff()) {
1316 Abort("InitPlaneEff", TSelector::kAbortProcess);
1320 if (strcmp(gProgName,"alieve") == 0)
1321 fRunAliEVE = InitAliEVE();
1326 //_____________________________________________________________________________
1327 Bool_t AliReconstruction::Process(Long64_t entry)
1329 // run the reconstruction over a single entry
1330 // from the chain with raw data
1331 AliCodeTimerAuto("");
1333 TTree *currTree = fChain->GetTree();
1334 AliRawEvent *event = new AliRawEvent;
1335 currTree->SetBranchAddress("rawevent",&event);
1336 currTree->GetEntry(entry);
1337 fRawReader = new AliRawReaderRoot(event);
1338 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1346 //_____________________________________________________________________________
1347 void AliReconstruction::Init(TTree *tree)
1350 AliError("The input tree is not found!");
1356 //_____________________________________________________________________________
1357 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1359 // run the reconstruction over a single event
1360 // The event loop is steered in Run method
1362 AliCodeTimerAuto("");
1364 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1365 fRunLoader->SetEventNumber(iEvent);
1366 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1368 fRunLoader->TreeE()->Fill();
1369 if (fRawReader && fRawReader->UseAutoSaveESD())
1370 fRunLoader->TreeE()->AutoSave("SaveSelf");
1373 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1377 AliInfo(Form("processing event %d", iEvent));
1379 // Fill Event-info object
1381 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1383 //Start of cycle for the in-loop QA
1384 if (fInLoopQA && fRunQA) {
1385 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam, fSameQACycle, fInLoopQA) ;
1387 if (fInLoopQA && fRunGlobalQA) {
1388 fSameQACycle = kFALSE;
1389 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1390 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1391 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1392 fSameQACycle = kTRUE;
1394 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1395 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1396 fSameQACycle = kTRUE;
1400 fRunLoader->GetEvent(iEvent);
1403 if (fInLoopQA && fRunQA)
1404 fQASteer->RunOneEvent(fRawReader) ;
1406 // local single event reconstruction
1407 if (!fRunLocalReconstruction.IsNull()) {
1408 TString detectors=fRunLocalReconstruction;
1409 // run HLT event reconstruction first
1410 // ;-( IsSelected changes the string
1411 if (IsSelected("HLT", detectors) &&
1412 !RunLocalEventReconstruction("HLT")) {
1413 if (fStopOnError) {CleanUp(); return kFALSE;}
1415 detectors=fRunLocalReconstruction;
1416 detectors.ReplaceAll("HLT", "");
1417 if (!RunLocalEventReconstruction(detectors)) {
1418 if (fStopOnError) {CleanUp(); return kFALSE;}
1422 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1423 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1424 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1425 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1427 // Set magnetic field from the tracker
1428 fesd->SetMagneticField(AliTracker::GetBz());
1429 fhltesd->SetMagneticField(AliTracker::GetBz());
1433 // Fill raw-data error log into the ESD
1434 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1437 if (fRunVertexFinder) {
1438 if (!RunVertexFinder(fesd)) {
1439 if (fStopOnError) {CleanUp(); return kFALSE;}
1444 if (!fRunTracking.IsNull()) {
1445 if (fRunMuonTracking) {
1446 if (!RunMuonTracking(fesd)) {
1447 if (fStopOnError) {CleanUp(); return kFALSE;}
1453 if (!fRunTracking.IsNull()) {
1454 if (!RunTracking(fesd)) {
1455 if (fStopOnError) {CleanUp(); return kFALSE;}
1460 if (!fFillESD.IsNull()) {
1461 TString detectors=fFillESD;
1462 // run HLT first and on hltesd
1463 // ;-( IsSelected changes the string
1464 if (IsSelected("HLT", detectors) &&
1465 !FillESD(fhltesd, "HLT")) {
1466 if (fStopOnError) {CleanUp(); return kFALSE;}
1469 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1470 if (detectors.Contains("ALL")) {
1472 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1473 detectors += fgkDetectorName[idet];
1477 detectors.ReplaceAll("HLT", "");
1478 if (!FillESD(fesd, detectors)) {
1479 if (fStopOnError) {CleanUp(); return kFALSE;}
1483 // fill Event header information from the RawEventHeader
1484 if (fRawReader){FillRawEventHeaderESD(fesd);}
1487 AliESDpid::MakePID(fesd);
1489 if (fFillTriggerESD) {
1490 if (!FillTriggerESD(fesd)) {
1491 if (fStopOnError) {CleanUp(); return kFALSE;}
1498 // Propagate track to the beam pipe (if not already done by ITS)
1500 const Int_t ntracks = fesd->GetNumberOfTracks();
1501 const Double_t kBz = fesd->GetMagneticField();
1502 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1505 UShort_t *selectedIdx=new UShort_t[ntracks];
1507 for (Int_t itrack=0; itrack<ntracks; itrack++){
1508 const Double_t kMaxStep = 5; //max step over the material
1511 AliESDtrack *track = fesd->GetTrack(itrack);
1512 if (!track) continue;
1514 AliExternalTrackParam *tpcTrack =
1515 (AliExternalTrackParam *)track->GetTPCInnerParam();
1519 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1522 Int_t n=trkArray.GetEntriesFast();
1523 selectedIdx[n]=track->GetID();
1524 trkArray.AddLast(tpcTrack);
1527 //Tracks refitted by ITS should already be at the SPD vertex
1528 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1531 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1532 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1537 // Improve the reconstructed primary vertex position using the tracks
1539 TObject *obj = fOptions.FindObject("ITS");
1541 TString optITS = obj->GetTitle();
1542 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1543 fRunVertexFinderTracks=kFALSE;
1545 if (fRunVertexFinderTracks) {
1546 // TPC + ITS primary vertex
1547 ftVertexer->SetITSrefitRequired();
1548 if(fDiamondProfile && fMeanVertexConstraint) {
1549 ftVertexer->SetVtxStart(fDiamondProfile);
1551 ftVertexer->SetConstraintOff();
1553 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1555 if (pvtx->GetStatus()) {
1556 fesd->SetPrimaryVertex(pvtx);
1557 for (Int_t i=0; i<ntracks; i++) {
1558 AliESDtrack *t = fesd->GetTrack(i);
1559 t->RelateToVertex(pvtx, kBz, kVeryBig);
1564 // TPC-only primary vertex
1565 ftVertexer->SetITSrefitNotRequired();
1566 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1567 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1569 ftVertexer->SetConstraintOff();
1571 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1573 if (pvtx->GetStatus()) {
1574 fesd->SetPrimaryVertexTPC(pvtx);
1575 for (Int_t i=0; i<ntracks; i++) {
1576 AliESDtrack *t = fesd->GetTrack(i);
1577 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1583 delete[] selectedIdx;
1585 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1590 AliV0vertexer vtxer;
1591 vtxer.Tracks2V0vertices(fesd);
1593 if (fRunCascadeFinder) {
1595 AliCascadeVertexer cvtxer;
1596 cvtxer.V0sTracks2CascadeVertices(fesd);
1601 if (fCleanESD) CleanESD(fesd);
1604 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1605 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1606 qadm->Exec(AliQA::kESDS, fesd);
1609 if (fWriteESDfriend) {
1610 fesdf->~AliESDfriend();
1611 new (fesdf) AliESDfriend(); // Reset...
1612 fesd->GetESDfriend(fesdf);
1616 // Auto-save the ESD tree in case of prompt reco @P2
1617 if (fRawReader && fRawReader->UseAutoSaveESD())
1618 ftree->AutoSave("SaveSelf");
1624 if (fRunAliEVE) RunAliEVE();
1628 if (fWriteESDfriend) {
1629 fesdf->~AliESDfriend();
1630 new (fesdf) AliESDfriend(); // Reset...
1633 ProcInfo_t ProcInfo;
1634 gSystem->GetProcInfo(&ProcInfo);
1635 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1638 // End of cycle for the in-loop
1639 if (fInLoopQA && fRunQA) {
1640 fQASteer->RunOneEvent(fesd) ;
1641 fQASteer->EndOfCycle() ;
1643 if (fInLoopQA && fRunGlobalQA) {
1644 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1646 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1647 qadm->EndOfCycle(AliQA::kRECPOINTS);
1648 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1649 qadm->EndOfCycle(AliQA::kESDS);
1655 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1656 if (fReconstructor[iDet])
1657 fReconstructor[iDet]->SetRecoParam(NULL);
1663 //_____________________________________________________________________________
1664 void AliReconstruction::SlaveTerminate()
1666 // Finalize the run on the slave side
1667 // Called after the exit
1668 // from the event loop
1669 AliCodeTimerAuto("");
1671 if (fIsNewRunLoader) { // galice.root didn't exist
1672 fRunLoader->WriteHeader("OVERWRITE");
1673 fRunLoader->CdGAFile();
1674 fRunLoader->Write(0, TObject::kOverwrite);
1677 ftree->GetUserInfo()->Add(fesd);
1678 fhlttree->GetUserInfo()->Add(fhltesd);
1680 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1681 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1683 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1684 cdbMapCopy->SetOwner(1);
1685 cdbMapCopy->SetName("cdbMap");
1686 TIter iter(cdbMap->GetTable());
1689 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1690 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1691 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1692 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1695 TList *cdbListCopy = new TList();
1696 cdbListCopy->SetOwner(1);
1697 cdbListCopy->SetName("cdbList");
1699 TIter iter2(cdbList);
1702 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1703 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1706 ftree->GetUserInfo()->Add(cdbMapCopy);
1707 ftree->GetUserInfo()->Add(cdbListCopy);
1710 if(fESDPar.Contains("ESD.par")){
1711 AliInfo("Attaching ESD.par to Tree");
1712 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1713 ftree->GetUserInfo()->Add(fn);
1719 if (fWriteESDfriend)
1720 ftree->SetBranchStatus("ESDfriend*",0);
1721 // we want to have only one tree version number
1722 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1725 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1726 if (fRunPlaneEff && !FinishPlaneEff()) {
1727 AliWarning("Finish PlaneEff evaluation failed");
1730 //Finish QA and end of cycle for out-of-loop QA
1731 if (!fInLoopQA && fRunQA)
1732 fQASteer->Run(fRunLocalReconstruction.Data(), AliQA::kNULLTASKINDEX, fSameQACycle) ;
1733 if (!fInLoopQA && fRunGlobalQA) {
1734 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1736 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1737 qadm->EndOfCycle(AliQA::kRECPOINTS);
1738 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1739 qadm->EndOfCycle(AliQA::kESDS);
1748 //_____________________________________________________________________________
1749 void AliReconstruction::Terminate()
1751 // Create tags for the events in the ESD tree (the ESD tree is always present)
1752 // In case of empty events the tags will contain dummy values
1753 AliCodeTimerAuto("");
1755 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1756 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1758 // Cleanup of CDB manager: cache and active storages!
1759 AliCDBManager::Instance()->ClearCache();
1762 //_____________________________________________________________________________
1763 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1765 // run the local reconstruction
1767 static Int_t eventNr=0;
1768 AliCodeTimerAuto("")
1770 TString detStr = detectors;
1771 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1772 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1773 AliReconstructor* reconstructor = GetReconstructor(iDet);
1774 if (!reconstructor) continue;
1775 AliLoader* loader = fLoader[iDet];
1776 // Matthias April 2008: temporary fix to run HLT reconstruction
1777 // although the HLT loader is missing
1778 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1780 reconstructor->Reconstruct(fRawReader, NULL);
1783 reconstructor->Reconstruct(dummy, NULL);
1788 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1791 // conversion of digits
1792 if (fRawReader && reconstructor->HasDigitConversion()) {
1793 AliInfo(Form("converting raw data digits into root objects for %s",
1794 fgkDetectorName[iDet]));
1795 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1796 // fgkDetectorName[iDet]));
1797 loader->LoadDigits("update");
1798 loader->CleanDigits();
1799 loader->MakeDigitsContainer();
1800 TTree* digitsTree = loader->TreeD();
1801 reconstructor->ConvertDigits(fRawReader, digitsTree);
1802 loader->WriteDigits("OVERWRITE");
1803 loader->UnloadDigits();
1805 // local reconstruction
1806 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1807 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1808 loader->LoadRecPoints("update");
1809 loader->CleanRecPoints();
1810 loader->MakeRecPointsContainer();
1811 TTree* clustersTree = loader->TreeR();
1812 if (fRawReader && !reconstructor->HasDigitConversion()) {
1813 reconstructor->Reconstruct(fRawReader, clustersTree);
1815 loader->LoadDigits("read");
1816 TTree* digitsTree = loader->TreeD();
1818 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1819 if (fStopOnError) return kFALSE;
1821 reconstructor->Reconstruct(digitsTree, clustersTree);
1823 loader->UnloadDigits();
1826 // In-loop QA for local reconstrucion
1827 TString detQAStr(fQADetectors) ;
1828 if (fRunQA && fInLoopQA)
1829 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1831 loader->WriteRecPoints("OVERWRITE");
1832 loader->UnloadRecPoints();
1833 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1835 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1836 AliError(Form("the following detectors were not found: %s",
1838 if (fStopOnError) return kFALSE;
1844 //_____________________________________________________________________________
1845 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1847 // run the barrel tracking
1849 AliCodeTimerAuto("")
1851 AliESDVertex* vertex = NULL;
1852 Double_t vtxPos[3] = {0, 0, 0};
1853 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1854 TArrayF mcVertex(3);
1855 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1856 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1857 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1861 AliInfo("running the ITS vertex finder");
1863 fLoader[0]->LoadRecPoints();
1864 TTree* cltree = fLoader[0]->TreeR();
1866 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1867 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1870 AliError("Can't get the ITS cluster tree");
1872 fLoader[0]->UnloadRecPoints();
1875 AliError("Can't get the ITS loader");
1878 AliWarning("Vertex not found");
1879 vertex = new AliESDVertex();
1880 vertex->SetName("default");
1883 vertex->SetName("reconstructed");
1887 AliInfo("getting the primary vertex from MC");
1888 vertex = new AliESDVertex(vtxPos, vtxErr);
1892 vertex->GetXYZ(vtxPos);
1893 vertex->GetSigmaXYZ(vtxErr);
1895 AliWarning("no vertex reconstructed");
1896 vertex = new AliESDVertex(vtxPos, vtxErr);
1898 esd->SetPrimaryVertexSPD(vertex);
1899 // if SPD multiplicity has been determined, it is stored in the ESD
1900 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1901 if(mult)esd->SetMultiplicity(mult);
1903 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1904 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1911 //_____________________________________________________________________________
1912 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1914 // run the HLT barrel tracking
1916 AliCodeTimerAuto("")
1919 AliError("Missing runLoader!");
1923 AliInfo("running HLT tracking");
1925 // Get a pointer to the HLT reconstructor
1926 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1927 if (!reconstructor) return kFALSE;
1930 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1931 TString detName = fgkDetectorName[iDet];
1932 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1933 reconstructor->SetOption(detName.Data());
1934 AliTracker *tracker = reconstructor->CreateTracker();
1936 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1937 if (fStopOnError) return kFALSE;
1941 Double_t vtxErr[3]={0.005,0.005,0.010};
1942 const AliESDVertex *vertex = esd->GetVertex();
1943 vertex->GetXYZ(vtxPos);
1944 tracker->SetVertex(vtxPos,vtxErr);
1946 fLoader[iDet]->LoadRecPoints("read");
1947 TTree* tree = fLoader[iDet]->TreeR();
1949 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1952 tracker->LoadClusters(tree);
1954 if (tracker->Clusters2Tracks(esd) != 0) {
1955 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1959 tracker->UnloadClusters();
1967 //_____________________________________________________________________________
1968 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1970 // run the muon spectrometer tracking
1972 AliCodeTimerAuto("")
1975 AliError("Missing runLoader!");
1978 Int_t iDet = 7; // for MUON
1980 AliInfo("is running...");
1982 // Get a pointer to the MUON reconstructor
1983 AliReconstructor *reconstructor = GetReconstructor(iDet);
1984 if (!reconstructor) return kFALSE;
1987 TString detName = fgkDetectorName[iDet];
1988 AliDebug(1, Form("%s tracking", detName.Data()));
1989 AliTracker *tracker = reconstructor->CreateTracker();
1991 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1996 fLoader[iDet]->LoadRecPoints("read");
1998 tracker->LoadClusters(fLoader[iDet]->TreeR());
2000 Int_t rv = tracker->Clusters2Tracks(esd);
2004 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2008 fLoader[iDet]->UnloadRecPoints();
2010 tracker->UnloadClusters();
2018 //_____________________________________________________________________________
2019 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2021 // run the barrel tracking
2022 static Int_t eventNr=0;
2023 AliCodeTimerAuto("")
2025 AliInfo("running tracking");
2027 //Fill the ESD with the T0 info (will be used by the TOF)
2028 if (fReconstructor[11] && fLoader[11]) {
2029 fLoader[11]->LoadRecPoints("READ");
2030 TTree *treeR = fLoader[11]->TreeR();
2031 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2034 // pass 1: TPC + ITS inwards
2035 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2036 if (!fTracker[iDet]) continue;
2037 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2040 fLoader[iDet]->LoadRecPoints("read");
2041 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2042 TTree* tree = fLoader[iDet]->TreeR();
2044 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2047 fTracker[iDet]->LoadClusters(tree);
2048 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2050 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2051 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2054 // preliminary PID in TPC needed by the ITS tracker
2056 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2057 AliESDpid::MakePID(esd);
2059 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2062 // pass 2: ALL backwards
2064 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2065 if (!fTracker[iDet]) continue;
2066 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2069 if (iDet > 1) { // all except ITS, TPC
2071 fLoader[iDet]->LoadRecPoints("read");
2072 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2073 tree = fLoader[iDet]->TreeR();
2075 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2078 fTracker[iDet]->LoadClusters(tree);
2079 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2083 if (iDet>1) // start filling residuals for the "outer" detectors
2084 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2086 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2087 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2092 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2093 fTracker[iDet]->UnloadClusters();
2094 fLoader[iDet]->UnloadRecPoints();
2096 // updated PID in TPC needed by the ITS tracker -MI
2098 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2099 AliESDpid::MakePID(esd);
2101 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2103 //stop filling residuals for the "outer" detectors
2104 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2106 // pass 3: TRD + TPC + ITS refit inwards
2108 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2109 if (!fTracker[iDet]) continue;
2110 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2113 if (iDet<2) // start filling residuals for TPC and ITS
2114 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2116 if (fTracker[iDet]->RefitInward(esd) != 0) {
2117 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2120 // run postprocessing
2121 if (fTracker[iDet]->PostProcess(esd) != 0) {
2122 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2125 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2128 // write space-points to the ESD in case alignment data output
2130 if (fWriteAlignmentData)
2131 WriteAlignmentData(esd);
2133 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2134 if (!fTracker[iDet]) continue;
2136 fTracker[iDet]->UnloadClusters();
2137 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2138 fLoader[iDet]->UnloadRecPoints();
2139 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2141 // stop filling residuals for TPC and ITS
2142 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2148 //_____________________________________________________________________________
2149 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2151 // Remove the data which are not needed for the physics analysis.
2154 Int_t nTracks=esd->GetNumberOfTracks();
2155 Int_t nV0s=esd->GetNumberOfV0s();
2157 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2159 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2160 Bool_t rc=esd->Clean(cleanPars);
2162 nTracks=esd->GetNumberOfTracks();
2163 nV0s=esd->GetNumberOfV0s();
2165 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2170 //_____________________________________________________________________________
2171 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2173 // fill the event summary data
2175 AliCodeTimerAuto("")
2176 static Int_t eventNr=0;
2177 TString detStr = detectors;
2179 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2180 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2181 AliReconstructor* reconstructor = GetReconstructor(iDet);
2182 if (!reconstructor) continue;
2183 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2184 TTree* clustersTree = NULL;
2185 if (fLoader[iDet]) {
2186 fLoader[iDet]->LoadRecPoints("read");
2187 clustersTree = fLoader[iDet]->TreeR();
2188 if (!clustersTree) {
2189 AliError(Form("Can't get the %s clusters tree",
2190 fgkDetectorName[iDet]));
2191 if (fStopOnError) return kFALSE;
2194 if (fRawReader && !reconstructor->HasDigitConversion()) {
2195 reconstructor->FillESD(fRawReader, clustersTree, esd);
2197 TTree* digitsTree = NULL;
2198 if (fLoader[iDet]) {
2199 fLoader[iDet]->LoadDigits("read");
2200 digitsTree = fLoader[iDet]->TreeD();
2202 AliError(Form("Can't get the %s digits tree",
2203 fgkDetectorName[iDet]));
2204 if (fStopOnError) return kFALSE;
2207 reconstructor->FillESD(digitsTree, clustersTree, esd);
2208 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2210 if (fLoader[iDet]) {
2211 fLoader[iDet]->UnloadRecPoints();
2215 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2216 AliError(Form("the following detectors were not found: %s",
2218 if (fStopOnError) return kFALSE;
2220 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2225 //_____________________________________________________________________________
2226 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2228 // Reads the trigger decision which is
2229 // stored in Trigger.root file and fills
2230 // the corresponding esd entries
2232 AliCodeTimerAuto("")
2234 AliInfo("Filling trigger information into the ESD");
2237 AliCTPRawStream input(fRawReader);
2238 if (!input.Next()) {
2239 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2242 if (esd->GetTriggerMask() != input.GetClassMask())
2243 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2244 input.GetClassMask(),esd->GetTriggerMask()));
2245 if (esd->GetOrbitNumber() != input.GetOrbitID())
2246 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2247 input.GetOrbitID(),esd->GetOrbitNumber()));
2248 if (esd->GetBunchCrossNumber() != input.GetBCID())
2249 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2250 input.GetBCID(),esd->GetBunchCrossNumber()));
2253 // Here one has to add the filling of trigger inputs and
2254 // interaction records
2264 //_____________________________________________________________________________
2265 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2268 // Filling information from RawReader Header
2271 if (!fRawReader) return kFALSE;
2273 AliInfo("Filling information from RawReader Header");
2275 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2276 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2277 esd->SetPeriodNumber(fRawReader->GetPeriod());
2279 esd->SetTimeStamp(fRawReader->GetTimestamp());
2280 esd->SetEventType(fRawReader->GetType());
2286 //_____________________________________________________________________________
2287 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2289 // check whether detName is contained in detectors
2290 // if yes, it is removed from detectors
2292 // check if all detectors are selected
2293 if ((detectors.CompareTo("ALL") == 0) ||
2294 detectors.BeginsWith("ALL ") ||
2295 detectors.EndsWith(" ALL") ||
2296 detectors.Contains(" ALL ")) {
2301 // search for the given detector
2302 Bool_t result = kFALSE;
2303 if ((detectors.CompareTo(detName) == 0) ||
2304 detectors.BeginsWith(detName+" ") ||
2305 detectors.EndsWith(" "+detName) ||
2306 detectors.Contains(" "+detName+" ")) {
2307 detectors.ReplaceAll(detName, "");
2311 // clean up the detectors string
2312 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2313 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2314 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2319 //_____________________________________________________________________________
2320 Bool_t AliReconstruction::InitRunLoader()
2322 // get or create the run loader
2324 if (gAlice) delete gAlice;
2327 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2328 // load all base libraries to get the loader classes
2329 TString libs = gSystem->GetLibraries();
2330 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2331 TString detName = fgkDetectorName[iDet];
2332 if (detName == "HLT") continue;
2333 if (libs.Contains("lib" + detName + "base.so")) continue;
2334 gSystem->Load("lib" + detName + "base.so");
2336 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2338 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2343 fRunLoader->CdGAFile();
2344 fRunLoader->LoadgAlice();
2346 //PH This is a temporary fix to give access to the kinematics
2347 //PH that is needed for the labels of ITS clusters
2348 fRunLoader->LoadHeader();
2349 fRunLoader->LoadKinematics();
2351 } else { // galice.root does not exist
2353 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2355 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2356 AliConfig::GetDefaultEventFolderName(),
2359 AliError(Form("could not create run loader in file %s",
2360 fGAliceFileName.Data()));
2364 fIsNewRunLoader = kTRUE;
2365 fRunLoader->MakeTree("E");
2367 if (fNumberOfEventsPerFile > 0)
2368 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2370 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2376 //_____________________________________________________________________________
2377 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2379 // get the reconstructor object and the loader for a detector
2381 if (fReconstructor[iDet]) {
2382 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2383 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2384 fReconstructor[iDet]->SetRecoParam(par);
2386 return fReconstructor[iDet];
2389 // load the reconstructor object
2390 TPluginManager* pluginManager = gROOT->GetPluginManager();
2391 TString detName = fgkDetectorName[iDet];
2392 TString recName = "Ali" + detName + "Reconstructor";
2394 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2396 AliReconstructor* reconstructor = NULL;
2397 // first check if a plugin is defined for the reconstructor
2398 TPluginHandler* pluginHandler =
2399 pluginManager->FindHandler("AliReconstructor", detName);
2400 // if not, add a plugin for it
2401 if (!pluginHandler) {
2402 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2403 TString libs = gSystem->GetLibraries();
2404 if (libs.Contains("lib" + detName + "base.so") ||
2405 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2406 pluginManager->AddHandler("AliReconstructor", detName,
2407 recName, detName + "rec", recName + "()");
2409 pluginManager->AddHandler("AliReconstructor", detName,
2410 recName, detName, recName + "()");
2412 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2414 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2415 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2417 if (reconstructor) {
2418 TObject* obj = fOptions.FindObject(detName.Data());
2419 if (obj) reconstructor->SetOption(obj->GetTitle());
2420 reconstructor->Init();
2421 fReconstructor[iDet] = reconstructor;
2424 // get or create the loader
2425 if (detName != "HLT") {
2426 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2427 if (!fLoader[iDet]) {
2428 AliConfig::Instance()
2429 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2431 // first check if a plugin is defined for the loader
2433 pluginManager->FindHandler("AliLoader", detName);
2434 // if not, add a plugin for it
2435 if (!pluginHandler) {
2436 TString loaderName = "Ali" + detName + "Loader";
2437 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2438 pluginManager->AddHandler("AliLoader", detName,
2439 loaderName, detName + "base",
2440 loaderName + "(const char*, TFolder*)");
2441 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2443 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2445 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2446 fRunLoader->GetEventFolder());
2448 if (!fLoader[iDet]) { // use default loader
2449 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2451 if (!fLoader[iDet]) {
2452 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2453 if (fStopOnError) return NULL;
2455 fRunLoader->AddLoader(fLoader[iDet]);
2456 fRunLoader->CdGAFile();
2457 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2458 fRunLoader->Write(0, TObject::kOverwrite);
2463 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2464 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2465 reconstructor->SetRecoParam(par);
2467 return reconstructor;
2470 //_____________________________________________________________________________
2471 Bool_t AliReconstruction::CreateVertexer()
2473 // create the vertexer
2476 AliReconstructor* itsReconstructor = GetReconstructor(0);
2477 if (itsReconstructor) {
2478 fVertexer = itsReconstructor->CreateVertexer();
2481 AliWarning("couldn't create a vertexer for ITS");
2482 if (fStopOnError) return kFALSE;
2488 //_____________________________________________________________________________
2489 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2491 // create the trackers
2493 TString detStr = detectors;
2494 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2495 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2496 AliReconstructor* reconstructor = GetReconstructor(iDet);
2497 if (!reconstructor) continue;
2498 TString detName = fgkDetectorName[iDet];
2499 if (detName == "HLT") {
2500 fRunHLTTracking = kTRUE;
2503 if (detName == "MUON") {
2504 fRunMuonTracking = kTRUE;
2509 fTracker[iDet] = reconstructor->CreateTracker();
2510 if (!fTracker[iDet] && (iDet < 7)) {
2511 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2512 if (fStopOnError) return kFALSE;
2514 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2520 //_____________________________________________________________________________
2521 void AliReconstruction::CleanUp()
2523 // delete trackers and the run loader and close and delete the file
2525 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2526 delete fReconstructor[iDet];
2527 fReconstructor[iDet] = NULL;
2528 fLoader[iDet] = NULL;
2529 delete fTracker[iDet];
2530 fTracker[iDet] = NULL;
2541 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2542 delete fDiamondProfile;
2543 fDiamondProfile = NULL;
2544 delete fDiamondProfileTPC;
2545 fDiamondProfileTPC = NULL;
2554 delete fParentRawReader;
2555 fParentRawReader=NULL;
2564 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2566 // Write space-points which are then used in the alignment procedures
2567 // For the moment only ITS, TPC, TRD and TOF
2569 Int_t ntracks = esd->GetNumberOfTracks();
2570 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2572 AliESDtrack *track = esd->GetTrack(itrack);
2575 for (Int_t iDet = 3; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2576 nsp += track->GetNcls(iDet);
2578 if (iDet==0) { // ITS "extra" clusters
2579 track->GetClusters(iDet,idx);
2580 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2585 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2586 track->SetTrackPointArray(sp);
2588 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2589 AliTracker *tracker = fTracker[iDet];
2590 if (!tracker) continue;
2591 Int_t nspdet = track->GetClusters(iDet,idx);
2593 if (iDet==0) // ITS "extra" clusters
2594 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2596 if (nspdet <= 0) continue;
2600 while (isp2 < nspdet) {
2601 Bool_t isvalid=kTRUE;
2603 Int_t index=idx[isp++];
2604 if (index < 0) continue;
2606 TString dets = fgkDetectorName[iDet];
2607 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2608 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2609 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2610 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2611 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2613 isvalid = tracker->GetTrackPoint(index,p);
2616 if (!isvalid) continue;
2617 sp->AddPoint(isptrack,&p); isptrack++;
2624 //_____________________________________________________________________________
2625 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2627 // The method reads the raw-data error log
2628 // accumulated within the rawReader.
2629 // It extracts the raw-data errors related to
2630 // the current event and stores them into
2631 // a TClonesArray inside the esd object.
2633 if (!fRawReader) return;
2635 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2637 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2639 if (iEvent != log->GetEventNumber()) continue;
2641 esd->AddRawDataErrorLog(log);
2646 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2647 // Dump a file content into a char in TNamed
2649 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2650 Int_t kBytes = (Int_t)in.tellg();
2651 printf("Size: %d \n",kBytes);
2654 char* memblock = new char [kBytes];
2655 in.seekg (0, ios::beg);
2656 in.read (memblock, kBytes);
2658 TString fData(memblock,kBytes);
2659 fn = new TNamed(pName,fData);
2660 printf("fData Size: %d \n",fData.Sizeof());
2661 printf("pName Size: %d \n",pName.Sizeof());
2662 printf("fn Size: %d \n",fn->Sizeof());
2666 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2672 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2673 // This is not really needed in AliReconstruction at the moment
2674 // but can serve as a template
2676 TList *fList = fTree->GetUserInfo();
2677 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2678 printf("fn Size: %d \n",fn->Sizeof());
2680 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2681 const char* cdata = fn->GetTitle();
2682 printf("fTmp Size %d\n",fTmp.Sizeof());
2684 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2685 printf("calculated size %d\n",size);
2686 ofstream out(pName.Data(),ios::out | ios::binary);
2687 out.write(cdata,size);
2692 //_____________________________________________________________________________
2693 void AliReconstruction::CheckQA()
2695 // check the QA of SIM for this run and remove the detectors
2696 // with status Fatal
2698 TString newRunLocalReconstruction ;
2699 TString newRunTracking ;
2700 TString newFillESD ;
2702 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2703 TString detName(AliQA::GetDetName(iDet)) ;
2704 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2705 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2706 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2708 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2709 fRunLocalReconstruction.Contains("ALL") ) {
2710 newRunLocalReconstruction += detName ;
2711 newRunLocalReconstruction += " " ;
2713 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2714 fRunTracking.Contains("ALL") ) {
2715 newRunTracking += detName ;
2716 newRunTracking += " " ;
2718 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2719 fFillESD.Contains("ALL") ) {
2720 newFillESD += detName ;
2725 fRunLocalReconstruction = newRunLocalReconstruction ;
2726 fRunTracking = newRunTracking ;
2727 fFillESD = newFillESD ;
2730 //_____________________________________________________________________________
2731 Int_t AliReconstruction::GetDetIndex(const char* detector)
2733 // return the detector index corresponding to detector
2735 for (index = 0; index < fgkNDetectors ; index++) {
2736 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2741 //_____________________________________________________________________________
2742 Bool_t AliReconstruction::FinishPlaneEff() {
2744 // Here execute all the necessary operationis, at the end of the tracking phase,
2745 // in case that evaluation of PlaneEfficiencies was required for some detector.
2746 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2748 // This Preliminary version works only FOR ITS !!!!!
2749 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2752 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2755 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2756 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2757 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2758 if(fTracker[iDet]) {
2759 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2760 TString name=planeeff->GetName();
2762 TFile* pefile = TFile::Open(name, "RECREATE");
2763 ret=(Bool_t)planeeff->Write();
2765 if(planeeff->GetCreateHistos()) {
2766 TString hname=planeeff->GetName();
2767 hname+="Histo.root";
2768 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2774 //_____________________________________________________________________________
2775 Bool_t AliReconstruction::InitPlaneEff() {
2777 // Here execute all the necessary operations, before of the tracking phase,
2778 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2779 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2780 // which should be updated/recalculated.
2782 // This Preliminary version will work only FOR ITS !!!!!
2783 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2786 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2788 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2792 //_____________________________________________________________________________
2793 Bool_t AliReconstruction::InitAliEVE()
2795 // This method should be called only in case
2796 // AliReconstruction is run
2797 // within the alieve environment.
2798 // It will initialize AliEVE in a way
2799 // so that it can visualize event processed
2800 // by AliReconstruction.
2801 // The return flag shows whenever the
2802 // AliEVE initialization was successful or not.
2805 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2806 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2807 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2809 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2810 gROOT->ProcessLine("alieve_online_init()");
2815 //_____________________________________________________________________________
2816 void AliReconstruction::RunAliEVE()
2818 // Runs AliEVE visualisation of
2819 // the current event.
2820 // Should be executed only after
2821 // successful initialization of AliEVE.
2823 AliInfo("Running AliEVE...");
2824 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2825 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2829 //_____________________________________________________________________________
2830 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2832 // Allows to run QA for a selected set of detectors
2833 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2834 // all selected detectors run the same selected tasks
2836 if (!detAndAction.Contains(":")) {
2837 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2841 Int_t colon = detAndAction.Index(":") ;
2842 fQADetectors = detAndAction(0, colon) ;
2843 if (fQADetectors.Contains("ALL") )
2844 fQADetectors = fFillESD ;
2845 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2846 if (fQATasks.Contains("ALL") ) {
2847 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2849 fQATasks.ToUpper() ;
2851 if ( fQATasks.Contains("RAW") )
2852 tempo = Form("%d ", AliQA::kRAWS) ;
2853 if ( fQATasks.Contains("RECPOINT") )
2854 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2855 if ( fQATasks.Contains("ESD") )
2856 tempo += Form("%d ", AliQA::kESDS) ;
2858 if (fQATasks.IsNull()) {
2859 AliInfo("No QA requested\n") ;
2864 TString tempo(fQATasks) ;
2865 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2866 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2867 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2868 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2873 //_____________________________________________________________________________
2874 Bool_t AliReconstruction::InitRecoParams()
2876 // The method accesses OCDB and retrieves all
2877 // the available reco-param objects from there.
2879 Bool_t isOK = kTRUE;
2881 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2883 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2884 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2888 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2890 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2891 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2893 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2897 TObject *recoParamObj = entry->GetObject();
2898 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2899 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2900 // Registering them in AliRecoParam
2901 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2903 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2904 // The detector has only onse set of reco parameters
2905 // Registering it in AliRecoParam
2906 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2907 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2908 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2911 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
2915 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
2924 //_____________________________________________________________________________
2925 Bool_t AliReconstruction::GetEventInfo()
2927 // Fill the event info object
2929 AliCodeTimerAuto("")
2931 AliCentralTrigger *aCTP = NULL;
2933 fEventInfo.SetEventType(fRawReader->GetType());
2935 ULong64_t mask = fRawReader->GetClassMask();
2936 fEventInfo.SetTriggerMask(mask);
2937 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
2938 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
2940 aCTP = new AliCentralTrigger();
2941 TString configstr("");
2942 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2943 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
2947 aCTP->SetClassMask(mask);
2948 aCTP->SetClusterMask(clmask);
2951 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
2953 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
2954 aCTP = fRunLoader->GetTrigger();
2955 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
2956 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
2959 AliWarning("No trigger can be loaded! The trigger information will not be used!");
2964 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2966 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
2967 if (fRawReader) delete aCTP;
2971 UChar_t clustmask = 0;
2973 ULong64_t trmask = fEventInfo.GetTriggerMask();
2974 const TObjArray& classesArray = config->GetClasses();
2975 Int_t nclasses = classesArray.GetEntriesFast();
2976 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
2977 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
2979 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2980 fesd->SetTriggerClass(trclass->GetName(),trindex);
2981 if (trmask & (1 << trindex)) {
2983 trclasses += trclass->GetName();
2985 clustmask |= trclass->GetCluster()->GetClusterMask();
2989 fEventInfo.SetTriggerClasses(trclasses);
2991 // Set the information in ESD
2992 fesd->SetTriggerMask(trmask);
2993 fesd->SetTriggerCluster(clustmask);
2995 if (!aCTP->CheckTriggeredDetectors()) {
2996 if (fRawReader) delete aCTP;
3000 if (fRawReader) delete aCTP;
3002 // We have to fill also the HLT decision here!!
3008 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3010 // Match the detector list found in the rec.C or the default 'ALL'
3011 // to the list found in the GRP (stored there by the shuttle PP which
3012 // gets the information from ECS)
3013 static TString resultList;
3014 TString detList = detectorList;
3018 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3019 if ((detectorMask >> iDet) & 0x1) {
3020 TString det = AliDAQ::OfflineModuleName(iDet);
3021 if ((detList.CompareTo("ALL") == 0) ||
3022 detList.BeginsWith("ALL ") ||
3023 detList.EndsWith(" ALL") ||
3024 detList.Contains(" ALL ") ||
3025 (detList.CompareTo(det) == 0) ||
3026 detList.BeginsWith(det) ||
3027 detList.EndsWith(det) ||
3028 detList.Contains( " "+det+" " )) {
3029 if (!resultList.EndsWith(det + " ")) {
3038 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3039 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3040 if ((detList.CompareTo("ALL") == 0) ||
3041 detList.BeginsWith("ALL ") ||
3042 detList.EndsWith(" ALL") ||
3043 detList.Contains(" ALL ") ||
3044 (detList.CompareTo(hltDet) == 0) ||
3045 detList.BeginsWith(hltDet) ||
3046 detList.EndsWith(hltDet) ||
3047 detList.Contains( " "+hltDet+" " )) {
3048 resultList += hltDet;
3052 return resultList.Data();
3056 //______________________________________________________________________________
3057 void AliReconstruction::Abort(const char *method, EAbort what)
3059 // Abort processing. If what = kAbortProcess, the Process() loop will be
3060 // aborted. If what = kAbortFile, the current file in a chain will be
3061 // aborted and the processing will continue with the next file, if there
3062 // is no next file then Process() will be aborted. Abort() can also be
3063 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3064 // the SlaveTerminate() and Terminate() are always called. The abort flag
3065 // can be checked in these methods using GetAbort().
3067 // The method is overwritten in AliReconstruction for better handling of
3068 // reco specific errors
3070 if (!fStopOnError) return;
3074 TString whyMess = method;
3075 whyMess += " failed! Aborting...";
3077 AliError(whyMess.Data());
3080 TString mess = "Abort";
3081 if (fAbort == kAbortProcess)
3082 mess = "AbortProcess";
3083 else if (fAbort == kAbortFile)
3086 Info(mess, whyMess.Data());