1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TPluginManager.h>
117 #include <TGeoManager.h>
118 #include <TLorentzVector.h>
121 #include <TObjArray.h>
125 #include <TProofOutputFile.h>
126 #include <TParameter.h>
128 #include "AliReconstruction.h"
129 #include "AliCodeTimer.h"
130 #include "AliReconstructor.h"
132 #include "AliRunLoader.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderDate.h"
136 #include "AliRawReaderRoot.h"
137 #include "AliRawEventHeaderBase.h"
138 #include "AliRawEvent.h"
139 #include "AliESDEvent.h"
140 #include "AliESDMuonTrack.h"
141 #include "AliESDfriend.h"
142 #include "AliESDVertex.h"
143 #include "AliESDcascade.h"
144 #include "AliESDkink.h"
145 #include "AliESDtrack.h"
146 #include "AliESDCaloCluster.h"
147 #include "AliESDCaloCells.h"
148 #include "AliMultiplicity.h"
149 #include "AliTracker.h"
150 #include "AliVertexer.h"
151 #include "AliVertexerTracks.h"
152 #include "AliV0vertexer.h"
153 #include "AliCascadeVertexer.h"
154 #include "AliHeader.h"
155 #include "AliGenEventHeader.h"
157 #include "AliESDpid.h"
158 #include "AliESDtrack.h"
159 #include "AliESDPmdTrack.h"
161 #include "AliESDTagCreator.h"
163 #include "AliGeomManager.h"
164 #include "AliTrackPointArray.h"
165 #include "AliCDBManager.h"
166 #include "AliCDBStorage.h"
167 #include "AliCDBEntry.h"
168 #include "AliAlignObj.h"
170 #include "AliCentralTrigger.h"
171 #include "AliTriggerConfiguration.h"
172 #include "AliTriggerClass.h"
173 #include "AliTriggerCluster.h"
174 #include "AliCTPRawStream.h"
176 #include "AliQADataMakerRec.h"
177 #include "AliGlobalQADataMaker.h"
179 #include "AliQADataMakerSteer.h"
181 #include "AliPlaneEff.h"
183 #include "AliSysInfo.h" // memory snapshots
184 #include "AliRawHLTManager.h"
186 #include "AliMagWrapCheb.h"
188 #include "AliDetectorRecoParam.h"
189 #include "AliGRPRecoParam.h"
190 #include "AliRunInfo.h"
191 #include "AliEventInfo.h"
195 ClassImp(AliReconstruction)
197 //_____________________________________________________________________________
198 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
200 //_____________________________________________________________________________
201 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
203 fUniformField(kFALSE),
204 fForcedFieldMap(NULL),
205 fRunVertexFinder(kTRUE),
206 fRunVertexFinderTracks(kTRUE),
207 fRunHLTTracking(kFALSE),
208 fRunMuonTracking(kFALSE),
210 fRunCascadeFinder(kTRUE),
211 fStopOnError(kFALSE),
212 fWriteAlignmentData(kFALSE),
213 fWriteESDfriend(kFALSE),
214 fFillTriggerESD(kTRUE),
222 fRunLocalReconstruction("ALL"),
226 fUseTrackingErrorsForAlignment(""),
227 fGAliceFileName(gAliceFilename),
232 fNumberOfEventsPerFile(1),
234 fLoadAlignFromCDB(kTRUE),
235 fLoadAlignData("ALL"),
243 fParentRawReader(NULL),
248 fDiamondProfile(NULL),
249 fDiamondProfileTPC(NULL),
250 fMeanVertexConstraint(kTRUE),
254 fAlignObjArray(NULL),
257 fInitCDBCalled(kFALSE),
258 fSetRunNumberFromDataCalled(kFALSE),
264 fSameQACycle(kFALSE),
266 fRunPlaneEff(kFALSE),
275 fIsNewRunLoader(kFALSE),
279 // create reconstruction object with default parameters
282 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
283 fReconstructor[iDet] = NULL;
284 fLoader[iDet] = NULL;
285 fTracker[iDet] = NULL;
287 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
288 fQACycles[iDet] = 999999 ;
289 fQAWriteExpert[iDet] = kFALSE ;
295 //_____________________________________________________________________________
296 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
298 fUniformField(rec.fUniformField),
299 fForcedFieldMap(NULL),
300 fRunVertexFinder(rec.fRunVertexFinder),
301 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
302 fRunHLTTracking(rec.fRunHLTTracking),
303 fRunMuonTracking(rec.fRunMuonTracking),
304 fRunV0Finder(rec.fRunV0Finder),
305 fRunCascadeFinder(rec.fRunCascadeFinder),
306 fStopOnError(rec.fStopOnError),
307 fWriteAlignmentData(rec.fWriteAlignmentData),
308 fWriteESDfriend(rec.fWriteESDfriend),
309 fFillTriggerESD(rec.fFillTriggerESD),
311 fCleanESD(rec.fCleanESD),
312 fV0DCAmax(rec.fV0DCAmax),
313 fV0CsPmin(rec.fV0CsPmin),
317 fRunLocalReconstruction(rec.fRunLocalReconstruction),
318 fRunTracking(rec.fRunTracking),
319 fFillESD(rec.fFillESD),
320 fLoadCDB(rec.fLoadCDB),
321 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
322 fGAliceFileName(rec.fGAliceFileName),
323 fRawInput(rec.fRawInput),
324 fEquipIdMap(rec.fEquipIdMap),
325 fFirstEvent(rec.fFirstEvent),
326 fLastEvent(rec.fLastEvent),
327 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
329 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
330 fLoadAlignData(rec.fLoadAlignData),
331 fESDPar(rec.fESDPar),
332 fUseHLTData(rec.fUseHLTData),
338 fParentRawReader(NULL),
340 fRecoParam(rec.fRecoParam),
343 fDiamondProfile(rec.fDiamondProfile),
344 fDiamondProfileTPC(rec.fDiamondProfileTPC),
345 fMeanVertexConstraint(rec.fMeanVertexConstraint),
349 fAlignObjArray(rec.fAlignObjArray),
350 fCDBUri(rec.fCDBUri),
352 fInitCDBCalled(rec.fInitCDBCalled),
353 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
354 fQADetectors(rec.fQADetectors),
356 fQATasks(rec.fQATasks),
358 fRunGlobalQA(rec.fRunGlobalQA),
359 fSameQACycle(rec.fSameQACycle),
360 fRunPlaneEff(rec.fRunPlaneEff),
369 fIsNewRunLoader(rec.fIsNewRunLoader),
375 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
376 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
378 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
379 fReconstructor[iDet] = NULL;
380 fLoader[iDet] = NULL;
381 fTracker[iDet] = NULL;
384 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
385 fQACycles[iDet] = rec.fQACycles[iDet];
386 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
389 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
390 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
394 //_____________________________________________________________________________
395 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
397 // assignment operator
398 // Used in PROOF mode
399 // Be very careful while modifing it!
400 // Simple rules to follow:
401 // for persistent data members - use their assignment operators
402 // for non-persistent ones - do nothing or take the default values from constructor
403 // TSelector members should not be touched
404 if(&rec == this) return *this;
406 fUniformField = rec.fUniformField;
407 fForcedFieldMap = NULL;
408 fRunVertexFinder = rec.fRunVertexFinder;
409 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
410 fRunHLTTracking = rec.fRunHLTTracking;
411 fRunMuonTracking = rec.fRunMuonTracking;
412 fRunV0Finder = rec.fRunV0Finder;
413 fRunCascadeFinder = rec.fRunCascadeFinder;
414 fStopOnError = rec.fStopOnError;
415 fWriteAlignmentData = rec.fWriteAlignmentData;
416 fWriteESDfriend = rec.fWriteESDfriend;
417 fFillTriggerESD = rec.fFillTriggerESD;
419 fCleanESD = rec.fCleanESD;
420 fV0DCAmax = rec.fV0DCAmax;
421 fV0CsPmin = rec.fV0CsPmin;
425 fRunLocalReconstruction = rec.fRunLocalReconstruction;
426 fRunTracking = rec.fRunTracking;
427 fFillESD = rec.fFillESD;
428 fLoadCDB = rec.fLoadCDB;
429 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
430 fGAliceFileName = rec.fGAliceFileName;
431 fRawInput = rec.fRawInput;
432 fEquipIdMap = rec.fEquipIdMap;
433 fFirstEvent = rec.fFirstEvent;
434 fLastEvent = rec.fLastEvent;
435 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
437 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
438 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
441 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
442 fLoadAlignData = rec.fLoadAlignData;
443 fESDPar = rec.fESDPar;
444 fUseHLTData = rec.fUseHLTData;
446 delete fRunInfo; fRunInfo = NULL;
447 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
449 fEventInfo = rec.fEventInfo;
453 fParentRawReader = NULL;
455 fRecoParam = rec.fRecoParam;
457 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
458 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
459 delete fLoader[iDet]; fLoader[iDet] = NULL;
460 delete fTracker[iDet]; fTracker[iDet] = NULL;
463 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
464 fQACycles[iDet] = rec.fQACycles[iDet];
465 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
469 delete fDiamondProfile; fDiamondProfile = NULL;
470 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
471 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
472 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
473 fMeanVertexConstraint = rec.fMeanVertexConstraint;
475 delete fGRPData; fGRPData = NULL;
476 if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
478 delete fAlignObjArray; fAlignObjArray = NULL;
481 fSpecCDBUri.Delete();
482 fInitCDBCalled = rec.fInitCDBCalled;
483 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
484 fQADetectors = rec.fQADetectors;
486 fQATasks = rec.fQATasks;
488 fRunGlobalQA = rec.fRunGlobalQA;
489 fSameQACycle = rec.fSameQACycle;
490 fRunPlaneEff = rec.fRunPlaneEff;
499 fIsNewRunLoader = rec.fIsNewRunLoader;
506 //_____________________________________________________________________________
507 AliReconstruction::~AliReconstruction()
512 delete fForcedFieldMap;
514 if (fAlignObjArray) {
515 fAlignObjArray->Delete();
516 delete fAlignObjArray;
518 fSpecCDBUri.Delete();
520 AliCodeTimer::Instance()->Print();
523 //_____________________________________________________________________________
524 void AliReconstruction::InitCDB()
526 // activate a default CDB storage
527 // First check if we have any CDB storage set, because it is used
528 // to retrieve the calibration and alignment constants
529 AliCodeTimerAuto("");
531 if (fInitCDBCalled) return;
532 fInitCDBCalled = kTRUE;
534 AliCDBManager* man = AliCDBManager::Instance();
535 if (man->IsDefaultStorageSet())
537 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
538 AliWarning("Default CDB storage has been already set !");
539 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
540 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
541 fCDBUri = man->GetDefaultStorage()->GetURI();
544 if (fCDBUri.Length() > 0)
546 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
547 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
548 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
550 fCDBUri="local://$ALICE_ROOT";
551 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
552 AliWarning("Default CDB storage not yet set !!!!");
553 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
554 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
557 man->SetDefaultStorage(fCDBUri);
560 // Now activate the detector specific CDB storage locations
561 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
562 TObject* obj = fSpecCDBUri[i];
564 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
565 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
566 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
567 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
569 AliSysInfo::AddStamp("InitCDB");
572 //_____________________________________________________________________________
573 void AliReconstruction::SetDefaultStorage(const char* uri) {
574 // Store the desired default CDB storage location
575 // Activate it later within the Run() method
581 //_____________________________________________________________________________
582 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
583 // Store a detector-specific CDB storage location
584 // Activate it later within the Run() method
586 AliCDBPath aPath(calibType);
587 if(!aPath.IsValid()){
588 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
589 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
590 if(!strcmp(calibType, fgkDetectorName[iDet])) {
591 aPath.SetPath(Form("%s/*", calibType));
592 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
596 if(!aPath.IsValid()){
597 AliError(Form("Not a valid path or detector: %s", calibType));
602 // // check that calibType refers to a "valid" detector name
603 // Bool_t isDetector = kFALSE;
604 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
605 // TString detName = fgkDetectorName[iDet];
606 // if(aPath.GetLevel0() == detName) {
607 // isDetector = kTRUE;
613 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
617 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
618 if (obj) fSpecCDBUri.Remove(obj);
619 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
623 //_____________________________________________________________________________
624 Bool_t AliReconstruction::SetRunNumberFromData()
626 // The method is called in Run() in order
627 // to set a correct run number.
628 // In case of raw data reconstruction the
629 // run number is taken from the raw data header
631 if (fSetRunNumberFromDataCalled) return kTRUE;
632 fSetRunNumberFromDataCalled = kTRUE;
634 AliCDBManager* man = AliCDBManager::Instance();
637 if(fRawReader->NextEvent()) {
638 if(man->GetRun() > 0) {
639 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
641 man->SetRun(fRawReader->GetRunNumber());
642 fRawReader->RewindEvents();
645 if(man->GetRun() > 0) {
646 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
649 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
655 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
657 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
662 // read run number from gAlice
663 if(rl->GetHeader()) {
664 man->SetRun(rl->GetHeader()->GetRun());
669 AliError("Neither run-loader header nor RawReader objects are found !");
681 //_____________________________________________________________________________
682 void AliReconstruction::SetCDBLock() {
683 // Set CDB lock: from now on it is forbidden to reset the run number
684 // or the default storage or to activate any further storage!
686 AliCDBManager::Instance()->SetLock(1);
689 //_____________________________________________________________________________
690 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
692 // Read the alignment objects from CDB.
693 // Each detector is supposed to have the
694 // alignment objects in DET/Align/Data CDB path.
695 // All the detector objects are then collected,
696 // sorted by geometry level (starting from ALIC) and
697 // then applied to the TGeo geometry.
698 // Finally an overlaps check is performed.
700 // Load alignment data from CDB and fill fAlignObjArray
701 if(fLoadAlignFromCDB){
703 TString detStr = detectors;
704 TString loadAlObjsListOfDets = "";
706 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
707 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
708 loadAlObjsListOfDets += fgkDetectorName[iDet];
709 loadAlObjsListOfDets += " ";
710 } // end loop over detectors
711 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
712 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
713 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
715 // Check if the array with alignment objects was
716 // provided by the user. If yes, apply the objects
717 // to the present TGeo geometry
718 if (fAlignObjArray) {
719 if (gGeoManager && gGeoManager->IsClosed()) {
720 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
721 AliError("The misalignment of one or more volumes failed!"
722 "Compare the list of simulated detectors and the list of detector alignment data!");
727 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
733 if (fAlignObjArray) {
734 fAlignObjArray->Delete();
735 delete fAlignObjArray; fAlignObjArray=NULL;
741 //_____________________________________________________________________________
742 void AliReconstruction::SetGAliceFile(const char* fileName)
744 // set the name of the galice file
746 fGAliceFileName = fileName;
749 //_____________________________________________________________________________
750 void AliReconstruction::SetInput(const char* input)
752 // In case the input string starts with 'mem://', we run in an online mode
753 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
754 // file is assumed. One can give as an input:
755 // mem://: - events taken from DAQ monitoring libs online
757 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
758 if (input) fRawInput = input;
761 //_____________________________________________________________________________
762 void AliReconstruction::SetOption(const char* detector, const char* option)
764 // set options for the reconstruction of a detector
766 TObject* obj = fOptions.FindObject(detector);
767 if (obj) fOptions.Remove(obj);
768 fOptions.Add(new TNamed(detector, option));
771 //_____________________________________________________________________________
772 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
774 // Set custom reconstruction parameters for a given detector
775 // Single set of parameters for all the events
777 // First check if the reco-params are global
778 if(!strcmp(detector, "GRP")) {
780 fRecoParam.AddDetRecoParam(fgkNDetectors,par);
784 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
785 if(!strcmp(detector, fgkDetectorName[iDet])) {
787 fRecoParam.AddDetRecoParam(iDet,par);
794 //_____________________________________________________________________________
795 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
796 //------------------------------------------------
797 // The magnetic field map, defined externally...
798 // L3 current 30000 A -> 0.5 T
799 // L3 current 12000 A -> 0.2 T
800 // dipole current 6000 A
801 // The polarities must be the same
802 //------------------------------------------------
803 const Float_t l3NominalCurrent1=30000.; // (A)
804 const Float_t l3NominalCurrent2=12000.; // (A)
805 const Float_t diNominalCurrent =6000. ; // (A)
807 const Float_t tolerance=0.03; // relative current tolerance
808 const Float_t zero=77.; // "zero" current (A)
811 Bool_t dipoleON=kFALSE;
813 TString s=(factor < 0) ? "L3: -" : "L3: +";
815 l3Current = TMath::Abs(l3Current);
816 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
817 map=AliMagWrapCheb::k5kG;
820 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
821 map=AliMagWrapCheb::k2kG;
824 if (l3Current < zero) {
825 map=AliMagWrapCheb::k2kG;
827 factor=0.; // in fact, this is a global factor...
828 fUniformField=kTRUE; // track with the uniform (zero) B field
830 AliError(Form("Wrong L3 current (%f A)!",l3Current));
834 diCurrent = TMath::Abs(diCurrent);
835 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
836 // 3% current tolerance...
840 if (diCurrent < zero) { // some small current..
844 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
848 delete fForcedFieldMap;
850 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
852 fForcedFieldMap->Print();
854 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
860 Bool_t AliReconstruction::InitGRP() {
861 //------------------------------------
862 // Initialization of the GRP entry
863 //------------------------------------
864 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
867 fGRPData = dynamic_cast<TMap*>(entry->GetObject());
869 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
873 AliError("No GRP entry found in OCDB!");
877 TObjString *lhcState=
878 dynamic_cast<TObjString*>(fGRPData->GetValue("fLHCState"));
880 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
883 TObjString *beamType=
884 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamType"));
886 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
889 TObjString *beamEnergyStr=
890 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamEnergy"));
891 if (!beamEnergyStr) {
892 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
896 dynamic_cast<TObjString*>(fGRPData->GetValue("fRunType"));
898 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
901 TObjString *activeDetectors=
902 dynamic_cast<TObjString*>(fGRPData->GetValue("fDetectorMask"));
903 if (!activeDetectors) {
904 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
907 fRunInfo = new AliRunInfo(lhcState ? lhcState->GetString().Data() : "UNKNOWN",
908 beamType ? beamType->GetString().Data() : "UNKNOWN",
909 beamEnergyStr ? beamEnergyStr->GetString().Atof() : 0,
910 runType ? runType->GetString().Data() : "UNKNOWN",
911 activeDetectors ? activeDetectors->GetString().Atoi() : 1074790399);
913 // Process the list of active detectors
914 if (activeDetectors && activeDetectors->GetString().IsDigit()) {
915 UInt_t detMask = activeDetectors->GetString().Atoi();
916 fLoadCDB.Form("%s %s %s %s",
917 fRunLocalReconstruction.Data(),
920 fQADetectors.Data());
921 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
922 fRunTracking = MatchDetectorList(fRunTracking,detMask);
923 fFillESD = MatchDetectorList(fFillESD,detMask);
924 fQADetectors = MatchDetectorList(fQADetectors,detMask);
925 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
928 AliInfo("===================================================================================");
929 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
930 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
931 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
932 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
933 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
934 AliInfo("===================================================================================");
936 //*** Dealing with the magnetic field map
937 if (AliTracker::GetFieldMap()) {
938 AliInfo("Running with the externally set B field !");
940 // Construct the field map out of the information retrieved from GRP.
945 TObjString *l3Current=
946 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
948 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
951 TObjString *l3Polarity=
952 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
954 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
959 TObjString *diCurrent=
960 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
962 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
965 TObjString *diPolarity=
966 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
968 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
973 Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
974 Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
975 Float_t l3Pol=atof(l3Polarity->GetName());
977 if (l3Pol != 0.) factor=-1.;
980 if (!SetFieldMap(l3Cur, diCur, factor)) {
981 AliFatal("Failed to creat a B field map ! Exiting...");
983 AliInfo("Running with the B field constructed out of GRP !");
986 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
991 //*** Get the diamond profile from OCDB
992 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
994 if (fMeanVertexConstraint)
995 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
997 AliError("No diamond profile found in OCDB!");
1000 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1002 if (fMeanVertexConstraint)
1003 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1005 AliError("No diamond profile found in OCDB!");
1011 //_____________________________________________________________________________
1012 Bool_t AliReconstruction::LoadCDB()
1014 AliCodeTimerAuto("");
1016 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1018 TString detStr = fLoadCDB;
1019 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1020 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1021 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1026 //_____________________________________________________________________________
1027 Bool_t AliReconstruction::Run(const char* input)
1030 AliCodeTimerAuto("");
1033 if (GetAbort() != TSelector::kContinue) return kFALSE;
1035 TChain *chain = NULL;
1036 if (fRawReader && (chain = fRawReader->GetChain())) {
1039 gProof->AddInput(this);
1041 outputFile.SetProtocol("root",kTRUE);
1042 outputFile.SetHost(gSystem->HostName());
1043 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1044 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1045 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1047 chain->Process("AliReconstruction");
1050 chain->Process(this);
1055 if (GetAbort() != TSelector::kContinue) return kFALSE;
1057 if (GetAbort() != TSelector::kContinue) return kFALSE;
1058 //******* The loop over events
1060 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1061 (fRawReader && fRawReader->NextEvent())) {
1062 if (!ProcessEvent(iEvent)) {
1063 Abort("ProcessEvent",TSelector::kAbortFile);
1069 if (GetAbort() != TSelector::kContinue) return kFALSE;
1071 if (GetAbort() != TSelector::kContinue) return kFALSE;
1077 //_____________________________________________________________________________
1078 void AliReconstruction::InitRawReader(const char* input)
1080 AliCodeTimerAuto("");
1082 // Init raw-reader and
1083 // set the input in case of raw data
1084 if (input) fRawInput = input;
1085 fRawReader = AliRawReader::Create(fRawInput.Data());
1087 AliInfo("Reconstruction will run over digits");
1089 if (!fEquipIdMap.IsNull() && fRawReader)
1090 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1092 if (!fUseHLTData.IsNull()) {
1093 // create the RawReaderHLT which performs redirection of HLT input data for
1094 // the specified detectors
1095 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1097 fParentRawReader=fRawReader;
1098 fRawReader=pRawReader;
1100 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1103 AliSysInfo::AddStamp("CreateRawReader");
1106 //_____________________________________________________________________________
1107 void AliReconstruction::InitRun(const char* input)
1109 // Initialization of raw-reader,
1110 // run number, CDB etc.
1111 AliCodeTimerAuto("");
1112 AliSysInfo::AddStamp("Start");
1114 // Initialize raw-reader if any
1115 InitRawReader(input);
1117 // Initialize the CDB storage
1120 // Set run number in CDBManager (if it is not already set by the user)
1121 if (!SetRunNumberFromData()) {
1122 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1126 // Set CDB lock: from now on it is forbidden to reset the run number
1127 // or the default storage or to activate any further storage!
1132 //_____________________________________________________________________________
1133 void AliReconstruction::Begin(TTree *)
1135 // Initialize AlReconstruction before
1136 // going into the event loop
1137 // Should follow the TSelector convention
1138 // i.e. initialize only the object on the client side
1139 AliCodeTimerAuto("");
1141 AliReconstruction *reco = NULL;
1143 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1146 AliSysInfo::AddStamp("ReadInputInBegin");
1149 // Import ideal TGeo geometry and apply misalignment
1151 TString geom(gSystem->DirName(fGAliceFileName));
1152 geom += "/geometry.root";
1153 AliGeomManager::LoadGeometry(geom.Data());
1155 Abort("LoadGeometry", TSelector::kAbortProcess);
1158 AliSysInfo::AddStamp("LoadGeom");
1159 TString detsToCheck=fRunLocalReconstruction;
1160 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1161 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1164 AliSysInfo::AddStamp("CheckGeom");
1167 if (!MisalignGeometry(fLoadAlignData)) {
1168 Abort("MisalignGeometry", TSelector::kAbortProcess);
1171 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1172 AliSysInfo::AddStamp("MisalignGeom");
1175 Abort("InitGRP", TSelector::kAbortProcess);
1178 AliSysInfo::AddStamp("InitGRP");
1181 Abort("LoadCDB", TSelector::kAbortProcess);
1184 AliSysInfo::AddStamp("LoadCDB");
1186 // Read the reconstruction parameters from OCDB
1187 if (!InitRecoParams()) {
1188 AliWarning("Not all detectors have correct RecoParam objects initialized");
1190 AliSysInfo::AddStamp("InitRecoParams");
1193 if (reco) *reco = *this;
1194 fInput->Add(gGeoManager);
1196 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1197 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1198 AliMagF *magFieldMap = (AliMagF*)AliTracker::GetFieldMap();
1199 magFieldMap->SetName("MagneticFieldMap");
1200 fInput->Add(magFieldMap);
1205 //_____________________________________________________________________________
1206 void AliReconstruction::SlaveBegin(TTree*)
1208 // Initialization related to run-loader,
1209 // vertexer, trackers, recontructors
1210 // In proof mode it is executed on the slave
1211 AliCodeTimerAuto("");
1213 TProofOutputFile *outProofFile = NULL;
1215 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1218 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1220 AliGeomManager::SetGeometry(tgeo);
1222 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1223 Int_t runNumber = -1;
1224 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1225 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1226 man->SetCacheFlag(kTRUE);
1227 man->SetLock(kTRUE);
1231 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1232 AliTracker::SetFieldMap(map,fUniformField);
1234 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1235 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1236 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1237 fOutput->Add(outProofFile);
1239 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1242 // get the run loader
1243 if (!InitRunLoader()) {
1244 Abort("InitRunLoader", TSelector::kAbortProcess);
1247 AliSysInfo::AddStamp("LoadLoader");
1249 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1250 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
1253 if (fRunVertexFinder && !CreateVertexer()) {
1254 Abort("CreateVertexer", TSelector::kAbortProcess);
1257 AliSysInfo::AddStamp("CreateVertexer");
1260 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1261 Abort("CreateTrackers", TSelector::kAbortProcess);
1264 AliSysInfo::AddStamp("CreateTrackers");
1266 // create the ESD output file and tree
1267 if (!outProofFile) {
1268 ffile = TFile::Open("AliESDs.root", "RECREATE");
1269 ffile->SetCompressionLevel(2);
1270 if (!ffile->IsOpen()) {
1271 Abort("OpenESDFile", TSelector::kAbortProcess);
1276 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1277 Abort(Form("Problems opening output PROOF file: %s/%s",
1278 outProofFile->GetDir(), outProofFile->GetFileName()),
1279 TSelector::kAbortProcess);
1284 ftree = new TTree("esdTree", "Tree with ESD objects");
1285 fesd = new AliESDEvent();
1286 fesd->CreateStdContent();
1287 fesd->WriteToTree(ftree);
1289 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1290 fhltesd = new AliESDEvent();
1291 fhltesd->CreateStdContent();
1292 fhltesd->WriteToTree(fhlttree);
1295 if (fWriteESDfriend) {
1296 fesdf = new AliESDfriend();
1297 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1298 br->SetFile("AliESDfriends.root");
1299 fesd->AddObject(fesdf);
1302 ProcInfo_t ProcInfo;
1303 gSystem->GetProcInfo(&ProcInfo);
1304 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1307 //Initialize the QA and start of cycle
1309 fQASteer = new AliQADataMakerSteer("rec") ;
1310 fQASteer->SetActiveDetectors(fQADetectors) ;
1311 for (Int_t det = 0 ; det < AliQA::kNDET ; det++) {
1312 fQASteer->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1313 if (fQAWriteExpert[det])
1314 fQASteer->SetWriteExpert(AliQA::DETECTORINDEX_t(det)) ;
1317 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1318 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1319 fQASteer->SetTasks(fQATasks) ;
1320 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam) ;
1324 Bool_t sameCycle = kFALSE ;
1325 if (!fQASteer) fQASteer = new AliQADataMakerSteer("rec") ;
1326 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1327 AliInfo(Form("Initializing the global QA data maker"));
1328 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1329 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1330 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS);
1331 AliTracker::SetResidualsArray(arr);
1334 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1335 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1336 qadm->Init(AliQA::kESDS);
1340 //Initialize the Plane Efficiency framework
1341 if (fRunPlaneEff && !InitPlaneEff()) {
1342 Abort("InitPlaneEff", TSelector::kAbortProcess);
1346 if (strcmp(gProgName,"alieve") == 0)
1347 fRunAliEVE = InitAliEVE();
1352 //_____________________________________________________________________________
1353 Bool_t AliReconstruction::Process(Long64_t entry)
1355 // run the reconstruction over a single entry
1356 // from the chain with raw data
1357 AliCodeTimerAuto("");
1359 TTree *currTree = fChain->GetTree();
1360 AliRawEvent *event = new AliRawEvent;
1361 currTree->SetBranchAddress("rawevent",&event);
1362 currTree->GetEntry(entry);
1363 fRawReader = new AliRawReaderRoot(event);
1364 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1372 //_____________________________________________________________________________
1373 void AliReconstruction::Init(TTree *tree)
1376 AliError("The input tree is not found!");
1382 //_____________________________________________________________________________
1383 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1385 // run the reconstruction over a single event
1386 // The event loop is steered in Run method
1388 AliCodeTimerAuto("");
1390 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1391 fRunLoader->SetEventNumber(iEvent);
1392 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1394 fRunLoader->TreeE()->Fill();
1395 if (fRawReader && fRawReader->UseAutoSaveESD())
1396 fRunLoader->TreeE()->AutoSave("SaveSelf");
1399 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1403 AliInfo(Form("processing event %d", iEvent));
1405 // Fill Event-info object
1407 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1409 // Set the reco-params
1411 TString detStr = fLoadCDB;
1412 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1413 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1414 AliReconstructor *reconstructor = GetReconstructor(iDet);
1415 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1416 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1417 reconstructor->SetRecoParam(par);
1422 fRunLoader->GetEvent(iEvent);
1426 fQASteer->RunOneEvent(fRawReader) ;
1428 // local single event reconstruction
1429 if (!fRunLocalReconstruction.IsNull()) {
1430 TString detectors=fRunLocalReconstruction;
1431 // run HLT event reconstruction first
1432 // ;-( IsSelected changes the string
1433 if (IsSelected("HLT", detectors) &&
1434 !RunLocalEventReconstruction("HLT")) {
1435 if (fStopOnError) {CleanUp(); return kFALSE;}
1437 detectors=fRunLocalReconstruction;
1438 detectors.ReplaceAll("HLT", "");
1439 if (!RunLocalEventReconstruction(detectors)) {
1440 if (fStopOnError) {CleanUp(); return kFALSE;}
1444 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1445 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1446 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1447 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1449 // Set magnetic field from the tracker
1450 fesd->SetMagneticField(AliTracker::GetBz());
1451 fhltesd->SetMagneticField(AliTracker::GetBz());
1455 // Fill raw-data error log into the ESD
1456 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1459 if (fRunVertexFinder) {
1460 if (!RunVertexFinder(fesd)) {
1461 if (fStopOnError) {CleanUp(); return kFALSE;}
1466 if (!fRunTracking.IsNull()) {
1467 if (fRunMuonTracking) {
1468 if (!RunMuonTracking(fesd)) {
1469 if (fStopOnError) {CleanUp(); return kFALSE;}
1475 if (!fRunTracking.IsNull()) {
1476 if (!RunTracking(fesd)) {
1477 if (fStopOnError) {CleanUp(); return kFALSE;}
1482 if (!fFillESD.IsNull()) {
1483 TString detectors=fFillESD;
1484 // run HLT first and on hltesd
1485 // ;-( IsSelected changes the string
1486 if (IsSelected("HLT", detectors) &&
1487 !FillESD(fhltesd, "HLT")) {
1488 if (fStopOnError) {CleanUp(); return kFALSE;}
1491 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1492 if (detectors.Contains("ALL")) {
1494 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1495 detectors += fgkDetectorName[idet];
1499 detectors.ReplaceAll("HLT", "");
1500 if (!FillESD(fesd, detectors)) {
1501 if (fStopOnError) {CleanUp(); return kFALSE;}
1505 // fill Event header information from the RawEventHeader
1506 if (fRawReader){FillRawEventHeaderESD(fesd);}
1509 AliESDpid::MakePID(fesd);
1511 if (fFillTriggerESD) {
1512 if (!FillTriggerESD(fesd)) {
1513 if (fStopOnError) {CleanUp(); return kFALSE;}
1520 // Propagate track to the beam pipe (if not already done by ITS)
1522 const Int_t ntracks = fesd->GetNumberOfTracks();
1523 const Double_t kBz = fesd->GetMagneticField();
1524 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1527 UShort_t *selectedIdx=new UShort_t[ntracks];
1529 for (Int_t itrack=0; itrack<ntracks; itrack++){
1530 const Double_t kMaxStep = 5; //max step over the material
1533 AliESDtrack *track = fesd->GetTrack(itrack);
1534 if (!track) continue;
1536 AliExternalTrackParam *tpcTrack =
1537 (AliExternalTrackParam *)track->GetTPCInnerParam();
1541 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1544 Int_t n=trkArray.GetEntriesFast();
1545 selectedIdx[n]=track->GetID();
1546 trkArray.AddLast(tpcTrack);
1549 //Tracks refitted by ITS should already be at the SPD vertex
1550 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1553 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1554 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1559 // Improve the reconstructed primary vertex position using the tracks
1561 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1562 if(fesd->GetPrimaryVertexSPD()) {
1563 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1564 if(vtitle.Contains("cosmics")) {
1565 runVertexFinderTracks=kFALSE;
1569 if (runVertexFinderTracks) {
1570 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1571 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(fgkNDetectors));
1573 // TPC + ITS primary vertex
1574 ftVertexer->SetITSMode();
1575 // get cuts for vertexer from AliGRPRecoParam
1577 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1578 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1579 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1580 ftVertexer->SetCuts(cutsVertexer);
1581 delete [] cutsVertexer; cutsVertexer = NULL;
1583 if(fDiamondProfile && fMeanVertexConstraint) {
1584 ftVertexer->SetVtxStart(fDiamondProfile);
1586 ftVertexer->SetConstraintOff();
1588 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1590 if (pvtx->GetStatus()) {
1591 fesd->SetPrimaryVertex(pvtx);
1592 for (Int_t i=0; i<ntracks; i++) {
1593 AliESDtrack *t = fesd->GetTrack(i);
1594 t->RelateToVertex(pvtx, kBz, kVeryBig);
1599 // TPC-only primary vertex
1600 ftVertexer->SetTPCMode();
1601 // get cuts for vertexer from AliGRPRecoParam
1603 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1604 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1605 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1606 ftVertexer->SetCuts(cutsVertexer);
1607 delete [] cutsVertexer; cutsVertexer = NULL;
1609 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1610 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1612 ftVertexer->SetConstraintOff();
1614 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1616 if (pvtx->GetStatus()) {
1617 fesd->SetPrimaryVertexTPC(pvtx);
1618 for (Int_t i=0; i<ntracks; i++) {
1619 AliESDtrack *t = fesd->GetTrack(i);
1620 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1626 delete[] selectedIdx;
1628 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1633 AliV0vertexer vtxer;
1634 vtxer.Tracks2V0vertices(fesd);
1636 if (fRunCascadeFinder) {
1638 AliCascadeVertexer cvtxer;
1639 cvtxer.V0sTracks2CascadeVertices(fesd);
1644 if (fCleanESD) CleanESD(fesd);
1647 fQASteer->RunOneEvent(fesd) ;
1650 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1651 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1652 qadm->Exec(AliQA::kESDS, fesd);
1655 if (fWriteESDfriend) {
1656 fesdf->~AliESDfriend();
1657 new (fesdf) AliESDfriend(); // Reset...
1658 fesd->GetESDfriend(fesdf);
1662 // Auto-save the ESD tree in case of prompt reco @P2
1663 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1664 ftree->AutoSave("SaveSelf");
1665 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1666 if (friendfile) friendfile->Save();
1673 if (fRunAliEVE) RunAliEVE();
1677 if (fWriteESDfriend) {
1678 fesdf->~AliESDfriend();
1679 new (fesdf) AliESDfriend(); // Reset...
1682 ProcInfo_t ProcInfo;
1683 gSystem->GetProcInfo(&ProcInfo);
1684 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1687 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1688 if (fReconstructor[iDet])
1689 fReconstructor[iDet]->SetRecoParam(NULL);
1692 fQASteer->Increment() ;
1696 //_____________________________________________________________________________
1697 void AliReconstruction::SlaveTerminate()
1699 // Finalize the run on the slave side
1700 // Called after the exit
1701 // from the event loop
1702 AliCodeTimerAuto("");
1704 if (fIsNewRunLoader) { // galice.root didn't exist
1705 fRunLoader->WriteHeader("OVERWRITE");
1706 fRunLoader->CdGAFile();
1707 fRunLoader->Write(0, TObject::kOverwrite);
1710 ftree->GetUserInfo()->Add(fesd);
1711 fhlttree->GetUserInfo()->Add(fhltesd);
1713 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1714 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1716 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1717 cdbMapCopy->SetOwner(1);
1718 cdbMapCopy->SetName("cdbMap");
1719 TIter iter(cdbMap->GetTable());
1722 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1723 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1724 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1725 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1728 TList *cdbListCopy = new TList();
1729 cdbListCopy->SetOwner(1);
1730 cdbListCopy->SetName("cdbList");
1732 TIter iter2(cdbList);
1735 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1736 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1739 ftree->GetUserInfo()->Add(cdbMapCopy);
1740 ftree->GetUserInfo()->Add(cdbListCopy);
1743 if(fESDPar.Contains("ESD.par")){
1744 AliInfo("Attaching ESD.par to Tree");
1745 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1746 ftree->GetUserInfo()->Add(fn);
1752 if (fWriteESDfriend)
1753 ftree->SetBranchStatus("ESDfriend*",0);
1754 // we want to have only one tree version number
1755 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1758 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1759 if (fRunPlaneEff && !FinishPlaneEff()) {
1760 AliWarning("Finish PlaneEff evaluation failed");
1763 // End of cycle for the in-loop
1765 fQASteer->EndOfCycle() ;
1768 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1770 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1771 qadm->EndOfCycle(AliQA::kRECPOINTS);
1772 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1773 qadm->EndOfCycle(AliQA::kESDS);
1781 //_____________________________________________________________________________
1782 void AliReconstruction::Terminate()
1784 // Create tags for the events in the ESD tree (the ESD tree is always present)
1785 // In case of empty events the tags will contain dummy values
1786 AliCodeTimerAuto("");
1788 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1789 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1791 // Cleanup of CDB manager: cache and active storages!
1792 AliCDBManager::Instance()->ClearCache();
1795 //_____________________________________________________________________________
1796 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1798 // run the local reconstruction
1800 static Int_t eventNr=0;
1801 AliCodeTimerAuto("")
1803 TString detStr = detectors;
1804 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1805 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1806 AliReconstructor* reconstructor = GetReconstructor(iDet);
1807 if (!reconstructor) continue;
1808 AliLoader* loader = fLoader[iDet];
1809 // Matthias April 2008: temporary fix to run HLT reconstruction
1810 // although the HLT loader is missing
1811 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1813 reconstructor->Reconstruct(fRawReader, NULL);
1816 reconstructor->Reconstruct(dummy, NULL);
1821 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1824 // conversion of digits
1825 if (fRawReader && reconstructor->HasDigitConversion()) {
1826 AliInfo(Form("converting raw data digits into root objects for %s",
1827 fgkDetectorName[iDet]));
1828 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1829 // fgkDetectorName[iDet]));
1830 loader->LoadDigits("update");
1831 loader->CleanDigits();
1832 loader->MakeDigitsContainer();
1833 TTree* digitsTree = loader->TreeD();
1834 reconstructor->ConvertDigits(fRawReader, digitsTree);
1835 loader->WriteDigits("OVERWRITE");
1836 loader->UnloadDigits();
1838 // local reconstruction
1839 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1840 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1841 loader->LoadRecPoints("update");
1842 loader->CleanRecPoints();
1843 loader->MakeRecPointsContainer();
1844 TTree* clustersTree = loader->TreeR();
1845 if (fRawReader && !reconstructor->HasDigitConversion()) {
1846 reconstructor->Reconstruct(fRawReader, clustersTree);
1848 loader->LoadDigits("read");
1849 TTree* digitsTree = loader->TreeD();
1851 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1852 if (fStopOnError) return kFALSE;
1854 reconstructor->Reconstruct(digitsTree, clustersTree);
1856 loader->UnloadDigits();
1859 TString detQAStr(fQADetectors) ;
1861 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1863 loader->WriteRecPoints("OVERWRITE");
1864 loader->UnloadRecPoints();
1865 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1867 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1868 AliError(Form("the following detectors were not found: %s",
1870 if (fStopOnError) return kFALSE;
1876 //_____________________________________________________________________________
1877 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1879 // run the barrel tracking
1881 AliCodeTimerAuto("")
1883 AliESDVertex* vertex = NULL;
1884 Double_t vtxPos[3] = {0, 0, 0};
1885 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1886 TArrayF mcVertex(3);
1887 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1888 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1889 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1893 AliInfo("running the ITS vertex finder");
1895 fLoader[0]->LoadRecPoints();
1896 TTree* cltree = fLoader[0]->TreeR();
1898 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1899 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1902 AliError("Can't get the ITS cluster tree");
1904 fLoader[0]->UnloadRecPoints();
1907 AliError("Can't get the ITS loader");
1910 AliWarning("Vertex not found");
1911 vertex = new AliESDVertex();
1912 vertex->SetName("default");
1915 vertex->SetName("reconstructed");
1919 AliInfo("getting the primary vertex from MC");
1920 vertex = new AliESDVertex(vtxPos, vtxErr);
1924 vertex->GetXYZ(vtxPos);
1925 vertex->GetSigmaXYZ(vtxErr);
1927 AliWarning("no vertex reconstructed");
1928 vertex = new AliESDVertex(vtxPos, vtxErr);
1930 esd->SetPrimaryVertexSPD(vertex);
1931 // if SPD multiplicity has been determined, it is stored in the ESD
1932 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1933 if(mult)esd->SetMultiplicity(mult);
1935 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1936 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1943 //_____________________________________________________________________________
1944 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1946 // run the HLT barrel tracking
1948 AliCodeTimerAuto("")
1951 AliError("Missing runLoader!");
1955 AliInfo("running HLT tracking");
1957 // Get a pointer to the HLT reconstructor
1958 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1959 if (!reconstructor) return kFALSE;
1962 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1963 TString detName = fgkDetectorName[iDet];
1964 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1965 reconstructor->SetOption(detName.Data());
1966 AliTracker *tracker = reconstructor->CreateTracker();
1968 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1969 if (fStopOnError) return kFALSE;
1973 Double_t vtxErr[3]={0.005,0.005,0.010};
1974 const AliESDVertex *vertex = esd->GetVertex();
1975 vertex->GetXYZ(vtxPos);
1976 tracker->SetVertex(vtxPos,vtxErr);
1978 fLoader[iDet]->LoadRecPoints("read");
1979 TTree* tree = fLoader[iDet]->TreeR();
1981 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1984 tracker->LoadClusters(tree);
1986 if (tracker->Clusters2Tracks(esd) != 0) {
1987 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1991 tracker->UnloadClusters();
1999 //_____________________________________________________________________________
2000 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2002 // run the muon spectrometer tracking
2004 AliCodeTimerAuto("")
2007 AliError("Missing runLoader!");
2010 Int_t iDet = 7; // for MUON
2012 AliInfo("is running...");
2014 // Get a pointer to the MUON reconstructor
2015 AliReconstructor *reconstructor = GetReconstructor(iDet);
2016 if (!reconstructor) return kFALSE;
2019 TString detName = fgkDetectorName[iDet];
2020 AliDebug(1, Form("%s tracking", detName.Data()));
2021 AliTracker *tracker = reconstructor->CreateTracker();
2023 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2028 fLoader[iDet]->LoadRecPoints("read");
2030 tracker->LoadClusters(fLoader[iDet]->TreeR());
2032 Int_t rv = tracker->Clusters2Tracks(esd);
2036 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2040 fLoader[iDet]->UnloadRecPoints();
2042 tracker->UnloadClusters();
2050 //_____________________________________________________________________________
2051 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2053 // run the barrel tracking
2054 static Int_t eventNr=0;
2055 AliCodeTimerAuto("")
2057 AliInfo("running tracking");
2059 //Fill the ESD with the T0 info (will be used by the TOF)
2060 if (fReconstructor[11] && fLoader[11]) {
2061 fLoader[11]->LoadRecPoints("READ");
2062 TTree *treeR = fLoader[11]->TreeR();
2063 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2066 // pass 1: TPC + ITS inwards
2067 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2068 if (!fTracker[iDet]) continue;
2069 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2072 fLoader[iDet]->LoadRecPoints("read");
2073 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2074 TTree* tree = fLoader[iDet]->TreeR();
2076 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2079 fTracker[iDet]->LoadClusters(tree);
2080 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2082 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2083 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2086 // preliminary PID in TPC needed by the ITS tracker
2088 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2089 AliESDpid::MakePID(esd);
2091 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2094 // pass 2: ALL backwards
2096 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2097 if (!fTracker[iDet]) continue;
2098 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2101 if (iDet > 1) { // all except ITS, TPC
2103 fLoader[iDet]->LoadRecPoints("read");
2104 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2105 tree = fLoader[iDet]->TreeR();
2107 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2110 fTracker[iDet]->LoadClusters(tree);
2111 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2115 if (iDet>1) // start filling residuals for the "outer" detectors
2116 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2118 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2119 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2124 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2125 fTracker[iDet]->UnloadClusters();
2126 fLoader[iDet]->UnloadRecPoints();
2128 // updated PID in TPC needed by the ITS tracker -MI
2130 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2131 AliESDpid::MakePID(esd);
2133 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2135 //stop filling residuals for the "outer" detectors
2136 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2138 // pass 3: TRD + TPC + ITS refit inwards
2140 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2141 if (!fTracker[iDet]) continue;
2142 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2145 if (iDet<2) // start filling residuals for TPC and ITS
2146 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2148 if (fTracker[iDet]->RefitInward(esd) != 0) {
2149 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2152 // run postprocessing
2153 if (fTracker[iDet]->PostProcess(esd) != 0) {
2154 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2157 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2160 // write space-points to the ESD in case alignment data output
2162 if (fWriteAlignmentData)
2163 WriteAlignmentData(esd);
2165 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2166 if (!fTracker[iDet]) continue;
2168 fTracker[iDet]->UnloadClusters();
2169 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2170 fLoader[iDet]->UnloadRecPoints();
2171 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2173 // stop filling residuals for TPC and ITS
2174 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2180 //_____________________________________________________________________________
2181 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2183 // Remove the data which are not needed for the physics analysis.
2186 Int_t nTracks=esd->GetNumberOfTracks();
2187 Int_t nV0s=esd->GetNumberOfV0s();
2189 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2191 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2192 Bool_t rc=esd->Clean(cleanPars);
2194 nTracks=esd->GetNumberOfTracks();
2195 nV0s=esd->GetNumberOfV0s();
2197 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2202 //_____________________________________________________________________________
2203 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2205 // fill the event summary data
2207 AliCodeTimerAuto("")
2208 static Int_t eventNr=0;
2209 TString detStr = detectors;
2211 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2212 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2213 AliReconstructor* reconstructor = GetReconstructor(iDet);
2214 if (!reconstructor) continue;
2215 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2216 TTree* clustersTree = NULL;
2217 if (fLoader[iDet]) {
2218 fLoader[iDet]->LoadRecPoints("read");
2219 clustersTree = fLoader[iDet]->TreeR();
2220 if (!clustersTree) {
2221 AliError(Form("Can't get the %s clusters tree",
2222 fgkDetectorName[iDet]));
2223 if (fStopOnError) return kFALSE;
2226 if (fRawReader && !reconstructor->HasDigitConversion()) {
2227 reconstructor->FillESD(fRawReader, clustersTree, esd);
2229 TTree* digitsTree = NULL;
2230 if (fLoader[iDet]) {
2231 fLoader[iDet]->LoadDigits("read");
2232 digitsTree = fLoader[iDet]->TreeD();
2234 AliError(Form("Can't get the %s digits tree",
2235 fgkDetectorName[iDet]));
2236 if (fStopOnError) return kFALSE;
2239 reconstructor->FillESD(digitsTree, clustersTree, esd);
2240 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2242 if (fLoader[iDet]) {
2243 fLoader[iDet]->UnloadRecPoints();
2247 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2248 AliError(Form("the following detectors were not found: %s",
2250 if (fStopOnError) return kFALSE;
2252 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2257 //_____________________________________________________________________________
2258 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2260 // Reads the trigger decision which is
2261 // stored in Trigger.root file and fills
2262 // the corresponding esd entries
2264 AliCodeTimerAuto("")
2266 AliInfo("Filling trigger information into the ESD");
2269 AliCTPRawStream input(fRawReader);
2270 if (!input.Next()) {
2271 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2274 if (esd->GetTriggerMask() != input.GetClassMask())
2275 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2276 input.GetClassMask(),esd->GetTriggerMask()));
2277 if (esd->GetOrbitNumber() != input.GetOrbitID())
2278 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2279 input.GetOrbitID(),esd->GetOrbitNumber()));
2280 if (esd->GetBunchCrossNumber() != input.GetBCID())
2281 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2282 input.GetBCID(),esd->GetBunchCrossNumber()));
2285 // Here one has to add the filling of trigger inputs and
2286 // interaction records
2296 //_____________________________________________________________________________
2297 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2300 // Filling information from RawReader Header
2303 if (!fRawReader) return kFALSE;
2305 AliInfo("Filling information from RawReader Header");
2307 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2308 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2309 esd->SetPeriodNumber(fRawReader->GetPeriod());
2311 esd->SetTimeStamp(fRawReader->GetTimestamp());
2312 esd->SetEventType(fRawReader->GetType());
2318 //_____________________________________________________________________________
2319 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2321 // check whether detName is contained in detectors
2322 // if yes, it is removed from detectors
2324 // check if all detectors are selected
2325 if ((detectors.CompareTo("ALL") == 0) ||
2326 detectors.BeginsWith("ALL ") ||
2327 detectors.EndsWith(" ALL") ||
2328 detectors.Contains(" ALL ")) {
2333 // search for the given detector
2334 Bool_t result = kFALSE;
2335 if ((detectors.CompareTo(detName) == 0) ||
2336 detectors.BeginsWith(detName+" ") ||
2337 detectors.EndsWith(" "+detName) ||
2338 detectors.Contains(" "+detName+" ")) {
2339 detectors.ReplaceAll(detName, "");
2343 // clean up the detectors string
2344 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2345 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2346 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2351 //_____________________________________________________________________________
2352 Bool_t AliReconstruction::InitRunLoader()
2354 // get or create the run loader
2356 if (gAlice) delete gAlice;
2359 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2360 // load all base libraries to get the loader classes
2361 TString libs = gSystem->GetLibraries();
2362 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2363 TString detName = fgkDetectorName[iDet];
2364 if (detName == "HLT") continue;
2365 if (libs.Contains("lib" + detName + "base.so")) continue;
2366 gSystem->Load("lib" + detName + "base.so");
2368 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2370 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2375 fRunLoader->CdGAFile();
2376 fRunLoader->LoadgAlice();
2378 //PH This is a temporary fix to give access to the kinematics
2379 //PH that is needed for the labels of ITS clusters
2380 fRunLoader->LoadHeader();
2381 fRunLoader->LoadKinematics();
2383 } else { // galice.root does not exist
2385 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2387 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2388 AliConfig::GetDefaultEventFolderName(),
2391 AliError(Form("could not create run loader in file %s",
2392 fGAliceFileName.Data()));
2396 fIsNewRunLoader = kTRUE;
2397 fRunLoader->MakeTree("E");
2399 if (fNumberOfEventsPerFile > 0)
2400 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2402 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2408 //_____________________________________________________________________________
2409 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2411 // get the reconstructor object and the loader for a detector
2413 if (fReconstructor[iDet]) {
2414 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2415 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2416 fReconstructor[iDet]->SetRecoParam(par);
2418 return fReconstructor[iDet];
2421 // load the reconstructor object
2422 TPluginManager* pluginManager = gROOT->GetPluginManager();
2423 TString detName = fgkDetectorName[iDet];
2424 TString recName = "Ali" + detName + "Reconstructor";
2426 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2428 AliReconstructor* reconstructor = NULL;
2429 // first check if a plugin is defined for the reconstructor
2430 TPluginHandler* pluginHandler =
2431 pluginManager->FindHandler("AliReconstructor", detName);
2432 // if not, add a plugin for it
2433 if (!pluginHandler) {
2434 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2435 TString libs = gSystem->GetLibraries();
2436 if (libs.Contains("lib" + detName + "base.so") ||
2437 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2438 pluginManager->AddHandler("AliReconstructor", detName,
2439 recName, detName + "rec", recName + "()");
2441 pluginManager->AddHandler("AliReconstructor", detName,
2442 recName, detName, recName + "()");
2444 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2446 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2447 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2449 if (reconstructor) {
2450 TObject* obj = fOptions.FindObject(detName.Data());
2451 if (obj) reconstructor->SetOption(obj->GetTitle());
2452 reconstructor->Init();
2453 fReconstructor[iDet] = reconstructor;
2456 // get or create the loader
2457 if (detName != "HLT") {
2458 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2459 if (!fLoader[iDet]) {
2460 AliConfig::Instance()
2461 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2463 // first check if a plugin is defined for the loader
2465 pluginManager->FindHandler("AliLoader", detName);
2466 // if not, add a plugin for it
2467 if (!pluginHandler) {
2468 TString loaderName = "Ali" + detName + "Loader";
2469 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2470 pluginManager->AddHandler("AliLoader", detName,
2471 loaderName, detName + "base",
2472 loaderName + "(const char*, TFolder*)");
2473 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2475 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2477 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2478 fRunLoader->GetEventFolder());
2480 if (!fLoader[iDet]) { // use default loader
2481 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2483 if (!fLoader[iDet]) {
2484 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2485 if (fStopOnError) return NULL;
2487 fRunLoader->AddLoader(fLoader[iDet]);
2488 fRunLoader->CdGAFile();
2489 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2490 fRunLoader->Write(0, TObject::kOverwrite);
2495 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2496 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2497 reconstructor->SetRecoParam(par);
2499 return reconstructor;
2502 //_____________________________________________________________________________
2503 Bool_t AliReconstruction::CreateVertexer()
2505 // create the vertexer
2508 AliReconstructor* itsReconstructor = GetReconstructor(0);
2509 if (itsReconstructor) {
2510 fVertexer = itsReconstructor->CreateVertexer();
2513 AliWarning("couldn't create a vertexer for ITS");
2514 if (fStopOnError) return kFALSE;
2520 //_____________________________________________________________________________
2521 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2523 // create the trackers
2525 TString detStr = detectors;
2526 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2527 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2528 AliReconstructor* reconstructor = GetReconstructor(iDet);
2529 if (!reconstructor) continue;
2530 TString detName = fgkDetectorName[iDet];
2531 if (detName == "HLT") {
2532 fRunHLTTracking = kTRUE;
2535 if (detName == "MUON") {
2536 fRunMuonTracking = kTRUE;
2541 fTracker[iDet] = reconstructor->CreateTracker();
2542 if (!fTracker[iDet] && (iDet < 7)) {
2543 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2544 if (fStopOnError) return kFALSE;
2546 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2552 //_____________________________________________________________________________
2553 void AliReconstruction::CleanUp()
2555 // delete trackers and the run loader and close and delete the file
2557 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2558 delete fReconstructor[iDet];
2559 fReconstructor[iDet] = NULL;
2560 fLoader[iDet] = NULL;
2561 delete fTracker[iDet];
2562 fTracker[iDet] = NULL;
2573 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2574 delete fDiamondProfile;
2575 fDiamondProfile = NULL;
2576 delete fDiamondProfileTPC;
2577 fDiamondProfileTPC = NULL;
2586 delete fParentRawReader;
2587 fParentRawReader=NULL;
2596 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2598 // Write space-points which are then used in the alignment procedures
2599 // For the moment only ITS, TPC, TRD and TOF
2601 Int_t ntracks = esd->GetNumberOfTracks();
2602 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2604 AliESDtrack *track = esd->GetTrack(itrack);
2607 for (Int_t iDet = 3; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2608 nsp += track->GetNcls(iDet);
2610 if (iDet==0) { // ITS "extra" clusters
2611 track->GetClusters(iDet,idx);
2612 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2617 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2618 track->SetTrackPointArray(sp);
2620 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2621 AliTracker *tracker = fTracker[iDet];
2622 if (!tracker) continue;
2623 Int_t nspdet = track->GetClusters(iDet,idx);
2625 if (iDet==0) // ITS "extra" clusters
2626 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2628 if (nspdet <= 0) continue;
2632 while (isp2 < nspdet) {
2633 Bool_t isvalid=kTRUE;
2635 Int_t index=idx[isp++];
2636 if (index < 0) continue;
2638 TString dets = fgkDetectorName[iDet];
2639 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2640 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2641 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2642 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2643 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2645 isvalid = tracker->GetTrackPoint(index,p);
2648 if (!isvalid) continue;
2649 sp->AddPoint(isptrack,&p); isptrack++;
2656 //_____________________________________________________________________________
2657 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2659 // The method reads the raw-data error log
2660 // accumulated within the rawReader.
2661 // It extracts the raw-data errors related to
2662 // the current event and stores them into
2663 // a TClonesArray inside the esd object.
2665 if (!fRawReader) return;
2667 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2669 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2671 if (iEvent != log->GetEventNumber()) continue;
2673 esd->AddRawDataErrorLog(log);
2678 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2679 // Dump a file content into a char in TNamed
2681 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2682 Int_t kBytes = (Int_t)in.tellg();
2683 printf("Size: %d \n",kBytes);
2686 char* memblock = new char [kBytes];
2687 in.seekg (0, ios::beg);
2688 in.read (memblock, kBytes);
2690 TString fData(memblock,kBytes);
2691 fn = new TNamed(pName,fData);
2692 printf("fData Size: %d \n",fData.Sizeof());
2693 printf("pName Size: %d \n",pName.Sizeof());
2694 printf("fn Size: %d \n",fn->Sizeof());
2698 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2704 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2705 // This is not really needed in AliReconstruction at the moment
2706 // but can serve as a template
2708 TList *fList = fTree->GetUserInfo();
2709 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2710 printf("fn Size: %d \n",fn->Sizeof());
2712 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2713 const char* cdata = fn->GetTitle();
2714 printf("fTmp Size %d\n",fTmp.Sizeof());
2716 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2717 printf("calculated size %d\n",size);
2718 ofstream out(pName.Data(),ios::out | ios::binary);
2719 out.write(cdata,size);
2724 //_____________________________________________________________________________
2725 void AliReconstruction::CheckQA()
2727 // check the QA of SIM for this run and remove the detectors
2728 // with status Fatal
2730 TString newRunLocalReconstruction ;
2731 TString newRunTracking ;
2732 TString newFillESD ;
2734 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2735 TString detName(AliQA::GetDetName(iDet)) ;
2736 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2737 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2738 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2740 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2741 fRunLocalReconstruction.Contains("ALL") ) {
2742 newRunLocalReconstruction += detName ;
2743 newRunLocalReconstruction += " " ;
2745 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2746 fRunTracking.Contains("ALL") ) {
2747 newRunTracking += detName ;
2748 newRunTracking += " " ;
2750 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2751 fFillESD.Contains("ALL") ) {
2752 newFillESD += detName ;
2757 fRunLocalReconstruction = newRunLocalReconstruction ;
2758 fRunTracking = newRunTracking ;
2759 fFillESD = newFillESD ;
2762 //_____________________________________________________________________________
2763 Int_t AliReconstruction::GetDetIndex(const char* detector)
2765 // return the detector index corresponding to detector
2767 for (index = 0; index < fgkNDetectors ; index++) {
2768 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2773 //_____________________________________________________________________________
2774 Bool_t AliReconstruction::FinishPlaneEff() {
2776 // Here execute all the necessary operationis, at the end of the tracking phase,
2777 // in case that evaluation of PlaneEfficiencies was required for some detector.
2778 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2780 // This Preliminary version works only FOR ITS !!!!!
2781 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2784 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2787 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2788 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2789 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2790 if(fTracker[iDet]) {
2791 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2792 TString name=planeeff->GetName();
2794 TFile* pefile = TFile::Open(name, "RECREATE");
2795 ret=(Bool_t)planeeff->Write();
2797 if(planeeff->GetCreateHistos()) {
2798 TString hname=planeeff->GetName();
2799 hname+="Histo.root";
2800 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2806 //_____________________________________________________________________________
2807 Bool_t AliReconstruction::InitPlaneEff() {
2809 // Here execute all the necessary operations, before of the tracking phase,
2810 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2811 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2812 // which should be updated/recalculated.
2814 // This Preliminary version will work only FOR ITS !!!!!
2815 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2818 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2820 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2824 //_____________________________________________________________________________
2825 Bool_t AliReconstruction::InitAliEVE()
2827 // This method should be called only in case
2828 // AliReconstruction is run
2829 // within the alieve environment.
2830 // It will initialize AliEVE in a way
2831 // so that it can visualize event processed
2832 // by AliReconstruction.
2833 // The return flag shows whenever the
2834 // AliEVE initialization was successful or not.
2837 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2838 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2839 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2841 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2842 gROOT->ProcessLine("alieve_online_init()");
2847 //_____________________________________________________________________________
2848 void AliReconstruction::RunAliEVE()
2850 // Runs AliEVE visualisation of
2851 // the current event.
2852 // Should be executed only after
2853 // successful initialization of AliEVE.
2855 AliInfo("Running AliEVE...");
2856 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2860 //_____________________________________________________________________________
2861 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2863 // Allows to run QA for a selected set of detectors
2864 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2865 // all selected detectors run the same selected tasks
2867 if (!detAndAction.Contains(":")) {
2868 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2872 Int_t colon = detAndAction.Index(":") ;
2873 fQADetectors = detAndAction(0, colon) ;
2874 if (fQADetectors.Contains("ALL") )
2875 fQADetectors = fFillESD ;
2876 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2877 if (fQATasks.Contains("ALL") ) {
2878 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2880 fQATasks.ToUpper() ;
2882 if ( fQATasks.Contains("RAW") )
2883 tempo = Form("%d ", AliQA::kRAWS) ;
2884 if ( fQATasks.Contains("RECPOINT") )
2885 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2886 if ( fQATasks.Contains("ESD") )
2887 tempo += Form("%d ", AliQA::kESDS) ;
2889 if (fQATasks.IsNull()) {
2890 AliInfo("No QA requested\n") ;
2895 TString tempo(fQATasks) ;
2896 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2897 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2898 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2899 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2904 //_____________________________________________________________________________
2905 Bool_t AliReconstruction::InitRecoParams()
2907 // The method accesses OCDB and retrieves all
2908 // the available reco-param objects from there.
2910 Bool_t isOK = kTRUE;
2912 TString detStr = fLoadCDB;
2913 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2915 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2917 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2918 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2922 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2924 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2925 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2927 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2931 TObject *recoParamObj = entry->GetObject();
2932 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2933 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2934 // Registering them in AliRecoParam
2935 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2937 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2938 // The detector has only onse set of reco parameters
2939 // Registering it in AliRecoParam
2940 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2941 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2942 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2945 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
2949 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
2953 if (AliDebugLevel() > 0) fRecoParam.Print();
2958 //_____________________________________________________________________________
2959 Bool_t AliReconstruction::GetEventInfo()
2961 // Fill the event info object
2963 AliCodeTimerAuto("")
2965 AliCentralTrigger *aCTP = NULL;
2967 fEventInfo.SetEventType(fRawReader->GetType());
2969 ULong64_t mask = fRawReader->GetClassMask();
2970 fEventInfo.SetTriggerMask(mask);
2971 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
2972 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
2974 aCTP = new AliCentralTrigger();
2975 TString configstr("");
2976 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2977 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
2981 aCTP->SetClassMask(mask);
2982 aCTP->SetClusterMask(clmask);
2985 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
2987 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
2988 aCTP = fRunLoader->GetTrigger();
2989 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
2990 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
2993 AliWarning("No trigger can be loaded! The trigger information will not be used!");
2998 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3000 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3001 if (fRawReader) delete aCTP;
3005 UChar_t clustmask = 0;
3007 ULong64_t trmask = fEventInfo.GetTriggerMask();
3008 const TObjArray& classesArray = config->GetClasses();
3009 Int_t nclasses = classesArray.GetEntriesFast();
3010 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3011 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3013 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
3014 fesd->SetTriggerClass(trclass->GetName(),trindex);
3015 if (trmask & (1 << trindex)) {
3017 trclasses += trclass->GetName();
3019 clustmask |= trclass->GetCluster()->GetClusterMask();
3023 fEventInfo.SetTriggerClasses(trclasses);
3025 // Set the information in ESD
3026 fesd->SetTriggerMask(trmask);
3027 fesd->SetTriggerCluster(clustmask);
3029 if (!aCTP->CheckTriggeredDetectors()) {
3030 if (fRawReader) delete aCTP;
3034 if (fRawReader) delete aCTP;
3036 // We have to fill also the HLT decision here!!
3042 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3044 // Match the detector list found in the rec.C or the default 'ALL'
3045 // to the list found in the GRP (stored there by the shuttle PP which
3046 // gets the information from ECS)
3047 static TString resultList;
3048 TString detList = detectorList;
3052 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3053 if ((detectorMask >> iDet) & 0x1) {
3054 TString det = AliDAQ::OfflineModuleName(iDet);
3055 if ((detList.CompareTo("ALL") == 0) ||
3056 detList.BeginsWith("ALL ") ||
3057 detList.EndsWith(" ALL") ||
3058 detList.Contains(" ALL ") ||
3059 (detList.CompareTo(det) == 0) ||
3060 detList.BeginsWith(det) ||
3061 detList.EndsWith(det) ||
3062 detList.Contains( " "+det+" " )) {
3063 if (!resultList.EndsWith(det + " ")) {
3072 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3073 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3074 if ((detList.CompareTo("ALL") == 0) ||
3075 detList.BeginsWith("ALL ") ||
3076 detList.EndsWith(" ALL") ||
3077 detList.Contains(" ALL ") ||
3078 (detList.CompareTo(hltDet) == 0) ||
3079 detList.BeginsWith(hltDet) ||
3080 detList.EndsWith(hltDet) ||
3081 detList.Contains( " "+hltDet+" " )) {
3082 resultList += hltDet;
3086 return resultList.Data();
3090 //______________________________________________________________________________
3091 void AliReconstruction::Abort(const char *method, EAbort what)
3093 // Abort processing. If what = kAbortProcess, the Process() loop will be
3094 // aborted. If what = kAbortFile, the current file in a chain will be
3095 // aborted and the processing will continue with the next file, if there
3096 // is no next file then Process() will be aborted. Abort() can also be
3097 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3098 // the SlaveTerminate() and Terminate() are always called. The abort flag
3099 // can be checked in these methods using GetAbort().
3101 // The method is overwritten in AliReconstruction for better handling of
3102 // reco specific errors
3104 if (!fStopOnError) return;
3108 TString whyMess = method;
3109 whyMess += " failed! Aborting...";
3111 AliError(whyMess.Data());
3114 TString mess = "Abort";
3115 if (fAbort == kAbortProcess)
3116 mess = "AbortProcess";
3117 else if (fAbort == kAbortFile)
3120 Info(mess, whyMess.Data());