]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.cxx
Removing variable shadowing
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // class for running the reconstruction                                      //
21 //                                                                           //
22 // Clusters and tracks are created for all detectors and all events by       //
23 // typing:                                                                   //
24 //                                                                           //
25 //   AliReconstruction rec;                                                  //
26 //   rec.Run();                                                              //
27 //                                                                           //
28 // The Run method returns kTRUE in case of successful execution.             //
29 //                                                                           //
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   //
32 //                                                                           //
33 //   rec.SetInput("...");                                                    //
34 //                                                                           //
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                               //
40 //                                                                           //
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                 //
44 //                                                                           //
45 //   rec.SetEventRange(..., ...);                                            //
46 //                                                                           //
47 // The index -1 (default) can be used for the last event to indicate no      //
48 // upper limit of the event range.                                           //
49 //                                                                           //
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):                     //
55 //                                                                           //
56 //   rec.SetNumberOfEventsPerFile(...);                                      //
57 //                                                                           //
58 //                                                                           //
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                                                         //
62 //                                                                           //
63 //   rec.SetGAliceFile("...");                                               //
64 //                                                                           //
65 // The local reconstruction can be switched on or off for individual         //
66 // detectors by                                                              //
67 //                                                                           //
68 //   rec.SetRunLocalReconstruction("...");                                   //
69 //                                                                           //
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.                                 //
73 //                                                                           //
74 // The reconstruction of the primary vertex position can be switched off by  //
75 //                                                                           //
76 //   rec.SetRunVertexFinder(kFALSE);                                         //
77 //                                                                           //
78 // The tracking and the creation of ESD tracks can be switched on for        //
79 // selected detectors by                                                     //
80 //                                                                           //
81 //   rec.SetRunTracking("...");                                              //
82 //                                                                           //
83 // Uniform/nonuniform field tracking switches (default: uniform field)       //
84 //                                                                           //
85 //   rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
86 //                                                                           //
87 // The filling of additional ESD information can be steered by               //
88 //                                                                           //
89 //   rec.SetFillESD("...");                                                  //
90 //                                                                           //
91 // Again, for both methods the string specifies the list of detectors.       //
92 // The default is "ALL".                                                     //
93 //                                                                           //
94 // The call of the shortcut method                                           //
95 //                                                                           //
96 //   rec.SetRunReconstruction("...");                                        //
97 //                                                                           //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and    //
99 // SetFillESD with the same detector selecting string as argument.           //
100 //                                                                           //
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.            //
103 //                                                                           //
104 // The input data of a detector can be replaced by the corresponding HLT     //
105 // data by calling (usual detector string)                                   //
106 // SetUseHLTData("...");                                                     //
107 //                                                                           //
108 // For debug purposes the method SetCheckPointLevel can be used. If the      //
109 // argument is greater than 0, files with ESD events will be written after   //
110 // selected steps of the reconstruction for each event:                      //
111 //   level 1: after tracking and after filling of ESD (final)                //
112 //   level 2: in addition after each tracking step                           //
113 //   level 3: in addition after the filling of ESD for each detector         //
114 // If a final check point file exists for an event, this event will be       //
115 // skipped in the reconstruction. The tracking and the filling of ESD for    //
116 // a detector will be skipped as well, if the corresponding check point      //
117 // file exists. The ESD event will then be loaded from the file instead.     //
118 //                                                                           //
119 ///////////////////////////////////////////////////////////////////////////////
120
121 #include <TArrayF.h>
122 #include <TFile.h>
123 #include <TList.h>
124 #include <TSystem.h>
125 #include <TROOT.h>
126 #include <TPluginManager.h>
127 #include <TGeoManager.h>
128 #include <TLorentzVector.h>
129 #include <TArrayS.h>
130 #include <TArrayD.h>
131 #include <TObjArray.h>
132
133 #include "AliReconstruction.h"
134 #include "AliCodeTimer.h"
135 #include "AliReconstructor.h"
136 #include "AliLog.h"
137 #include "AliRunLoader.h"
138 #include "AliRun.h"
139 #include "AliRawReaderFile.h"
140 #include "AliRawReaderDate.h"
141 #include "AliRawReaderRoot.h"
142 #include "AliRawEventHeaderBase.h"
143 #include "AliESDEvent.h"
144 #include "AliESDMuonTrack.h"
145 #include "AliESDfriend.h"
146 #include "AliESDVertex.h"
147 #include "AliESDcascade.h"
148 #include "AliESDkink.h"
149 #include "AliESDtrack.h"
150 #include "AliESDCaloCluster.h"
151 #include "AliESDCaloCells.h"
152 #include "AliMultiplicity.h"
153 #include "AliTracker.h"
154 #include "AliVertexer.h"
155 #include "AliVertexerTracks.h"
156 #include "AliV0vertexer.h"
157 #include "AliCascadeVertexer.h"
158 #include "AliHeader.h"
159 #include "AliGenEventHeader.h"
160 #include "AliPID.h"
161 #include "AliESDpid.h"
162 #include "AliESDtrack.h"
163 #include "AliESDPmdTrack.h"
164
165 #include "AliESDTagCreator.h"
166 #include "AliAODTagCreator.h"
167
168 #include "AliGeomManager.h"
169 #include "AliTrackPointArray.h"
170 #include "AliCDBManager.h"
171 #include "AliCDBStorage.h"
172 #include "AliCDBEntry.h"
173 #include "AliAlignObj.h"
174
175 #include "AliCentralTrigger.h"
176 #include "AliTriggerConfiguration.h"
177 #include "AliTriggerClass.h"
178 #include "AliCTPRawStream.h"
179
180 #include "AliQADataMakerRec.h" 
181 #include "AliGlobalQADataMaker.h" 
182 #include "AliQA.h"
183 #include "AliQADataMakerSteer.h"
184
185 #include "AliPlaneEff.h"
186
187 #include "AliSysInfo.h" // memory snapshots
188 #include "AliRawHLTManager.h"
189
190
191 ClassImp(AliReconstruction)
192
193
194 //_____________________________________________________________________________
195 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
196
197 //_____________________________________________________________________________
198 AliReconstruction::AliReconstruction(const char* gAliceFilename,
199                                      const char* name, const char* title) :
200   TNamed(name, title),
201
202   fUniformField(kTRUE),
203   fRunVertexFinder(kTRUE),
204   fRunVertexFinderTracks(kTRUE),
205   fRunHLTTracking(kFALSE),
206   fRunMuonTracking(kFALSE),
207   fRunV0Finder(kTRUE),
208   fRunCascadeFinder(kTRUE),
209   fStopOnError(kFALSE),
210   fWriteAlignmentData(kFALSE),
211   fWriteESDfriend(kFALSE),
212   fWriteAOD(kFALSE),
213   fFillTriggerESD(kTRUE),
214
215   fCleanESD(kTRUE),
216   fV0DCAmax(3.),
217   fV0CsPmin(0.),
218   fDmax(50.),
219   fZmax(50.),
220
221   fRunLocalReconstruction("ALL"),
222   fRunTracking("ALL"),
223   fFillESD("ALL"),
224   fUseTrackingErrorsForAlignment(""),
225   fGAliceFileName(gAliceFilename),
226   fInput(""),
227   fEquipIdMap(""),
228   fFirstEvent(0),
229   fLastEvent(-1),
230   fNumberOfEventsPerFile(1),
231   fCheckPointLevel(0),
232   fOptions(),
233   fLoadAlignFromCDB(kTRUE),
234   fLoadAlignData("ALL"),
235   fESDPar(""),
236   fUseHLTData(),
237
238   fRunLoader(NULL),
239   fRawReader(NULL),
240   fParentRawReader(NULL),
241
242   fVertexer(NULL),
243   fDiamondProfile(NULL),
244   fDiamondProfileTPC(NULL),
245   fMeanVertexConstraint(kTRUE),
246
247   fGRPData(NULL),
248
249   fAlignObjArray(NULL),
250   fCDBUri(),
251   fSpecCDBUri(), 
252   fInitCDBCalled(kFALSE),
253   fSetRunNumberFromDataCalled(kFALSE),
254   fRunQA(kTRUE),  
255   fRunGlobalQA(kTRUE),
256   fInLoopQA(kFALSE),
257   fSameQACycle(kFALSE),
258
259   fRunPlaneEff(kFALSE),
260
261   fesd(NULL),
262   fhltesd(NULL),
263   fesdf(NULL),
264   ffile(NULL),
265   ftree(NULL),
266   fhlttree(NULL),
267   ffileOld(NULL),
268   ftreeOld(NULL),
269   fhlttreeOld(NULL),
270   ftVertexer(NULL),
271   fIsNewRunLoader(kFALSE),
272   fRunAliEVE(kFALSE)
273 {
274 // create reconstruction object with default parameters
275   
276   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
277     fReconstructor[iDet] = NULL;
278     fLoader[iDet] = NULL;
279     fTracker[iDet] = NULL;
280     fQADataMaker[iDet] = NULL;
281         fQACycles[iDet] = 999999;       
282   }
283   fQADataMaker[fgkNDetectors]=NULL;  //Global QA
284   AliPID pid;
285 }
286
287 //_____________________________________________________________________________
288 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
289   TNamed(rec),
290
291   fUniformField(rec.fUniformField),
292   fRunVertexFinder(rec.fRunVertexFinder),
293   fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
294   fRunHLTTracking(rec.fRunHLTTracking),
295   fRunMuonTracking(rec.fRunMuonTracking),
296   fRunV0Finder(rec.fRunV0Finder),
297   fRunCascadeFinder(rec.fRunCascadeFinder),
298   fStopOnError(rec.fStopOnError),
299   fWriteAlignmentData(rec.fWriteAlignmentData),
300   fWriteESDfriend(rec.fWriteESDfriend),
301   fWriteAOD(rec.fWriteAOD),
302   fFillTriggerESD(rec.fFillTriggerESD),
303
304   fCleanESD(rec.fCleanESD),
305   fV0DCAmax(rec.fV0DCAmax),
306   fV0CsPmin(rec.fV0CsPmin),
307   fDmax(rec.fDmax),
308   fZmax(rec.fZmax),
309
310   fRunLocalReconstruction(rec.fRunLocalReconstruction),
311   fRunTracking(rec.fRunTracking),
312   fFillESD(rec.fFillESD),
313   fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
314   fGAliceFileName(rec.fGAliceFileName),
315   fInput(rec.fInput),
316   fEquipIdMap(rec.fEquipIdMap),
317   fFirstEvent(rec.fFirstEvent),
318   fLastEvent(rec.fLastEvent),
319   fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
320   fCheckPointLevel(0),
321   fOptions(),
322   fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
323   fLoadAlignData(rec.fLoadAlignData),
324   fESDPar(rec.fESDPar),
325   fUseHLTData(rec.fUseHLTData),
326
327   fRunLoader(NULL),
328   fRawReader(NULL),
329   fParentRawReader(NULL),
330
331   fVertexer(NULL),
332   fDiamondProfile(NULL),
333   fDiamondProfileTPC(NULL),
334   fMeanVertexConstraint(rec.fMeanVertexConstraint),
335
336   fGRPData(NULL),
337
338   fAlignObjArray(rec.fAlignObjArray),
339   fCDBUri(rec.fCDBUri),
340   fSpecCDBUri(), 
341   fInitCDBCalled(rec.fInitCDBCalled),
342   fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
343   fRunQA(rec.fRunQA),  
344   fRunGlobalQA(rec.fRunGlobalQA),
345   fInLoopQA(rec.fInLoopQA),
346   fSameQACycle(rec.fSameQACycle),
347   fRunPlaneEff(rec.fRunPlaneEff),
348
349   fesd(NULL),
350   fhltesd(NULL),
351   fesdf(NULL),
352   ffile(NULL),
353   ftree(NULL),
354   fhlttree(NULL),
355   ffileOld(NULL),
356   ftreeOld(NULL),
357   fhlttreeOld(NULL),
358   ftVertexer(NULL),
359   fIsNewRunLoader(rec.fIsNewRunLoader),
360   fRunAliEVE(kFALSE)
361 {
362 // copy constructor
363
364   for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
365     if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
366   }
367   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
368     fReconstructor[iDet] = NULL;
369     fLoader[iDet] = NULL;
370     fTracker[iDet] = NULL;
371     fQADataMaker[iDet] = NULL;
372         fQACycles[iDet] = rec.fQACycles[iDet];  
373   }
374   fQADataMaker[fgkNDetectors]=NULL;  //Global QA
375   for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
376     if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
377   }
378 }
379
380 //_____________________________________________________________________________
381 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
382 {
383 // assignment operator
384
385   this->~AliReconstruction();
386   new(this) AliReconstruction(rec);
387   return *this;
388 }
389
390 //_____________________________________________________________________________
391 AliReconstruction::~AliReconstruction()
392 {
393 // clean up
394
395   CleanUp();
396   fOptions.Delete();
397   fSpecCDBUri.Delete();
398
399   AliCodeTimer::Instance()->Print();
400 }
401
402 //_____________________________________________________________________________
403 void AliReconstruction::InitCDB()
404 {
405 // activate a default CDB storage
406 // First check if we have any CDB storage set, because it is used 
407 // to retrieve the calibration and alignment constants
408
409   if (fInitCDBCalled) return;
410   fInitCDBCalled = kTRUE;
411
412   AliCDBManager* man = AliCDBManager::Instance();
413   if (man->IsDefaultStorageSet())
414   {
415     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
416     AliWarning("Default CDB storage has been already set !");
417     AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
418     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
419     fCDBUri = man->GetDefaultStorage()->GetURI();
420   }
421   else {
422     if (fCDBUri.Length() > 0) 
423     {
424         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
425         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
426         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
427     } else {
428         fCDBUri="local://$ALICE_ROOT";
429         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
430         AliWarning("Default CDB storage not yet set !!!!");
431         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
432         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
433                 
434     }
435     man->SetDefaultStorage(fCDBUri);
436   }
437
438   // Now activate the detector specific CDB storage locations
439   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
440     TObject* obj = fSpecCDBUri[i];
441     if (!obj) continue;
442     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
443     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
444     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
445     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
446   }
447   
448 }
449
450 //_____________________________________________________________________________
451 void AliReconstruction::SetDefaultStorage(const char* uri) {
452 // Store the desired default CDB storage location
453 // Activate it later within the Run() method
454
455   fCDBUri = uri;
456
457 }
458
459 //_____________________________________________________________________________
460 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
461 // Store a detector-specific CDB storage location
462 // Activate it later within the Run() method
463
464   AliCDBPath aPath(calibType);
465   if(!aPath.IsValid()){
466         // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
467         for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
468                 if(!strcmp(calibType, fgkDetectorName[iDet])) {
469                         aPath.SetPath(Form("%s/*", calibType));
470                         AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
471                         break;
472                 }
473         }
474         if(!aPath.IsValid()){
475                 AliError(Form("Not a valid path or detector: %s", calibType));
476                 return;
477         }
478   }
479
480 //  // check that calibType refers to a "valid" detector name
481 //  Bool_t isDetector = kFALSE;
482 //  for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
483 //    TString detName = fgkDetectorName[iDet];
484 //    if(aPath.GetLevel0() == detName) {
485 //      isDetector = kTRUE;
486 //      break;
487 //    }
488 //  }
489 //
490 //  if(!isDetector) {
491 //      AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
492 //      return;
493 //  }
494
495   TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
496   if (obj) fSpecCDBUri.Remove(obj);
497   fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
498
499 }
500
501 //_____________________________________________________________________________
502 Bool_t AliReconstruction::SetRunNumberFromData()
503 {
504   // The method is called in Run() in order
505   // to set a correct run number.
506   // In case of raw data reconstruction the
507   // run number is taken from the raw data header
508
509   if (fSetRunNumberFromDataCalled) return kTRUE;
510   fSetRunNumberFromDataCalled = kTRUE;
511   
512   AliCDBManager* man = AliCDBManager::Instance();
513   
514   if(man->GetRun() > 0) {
515         AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
516   } 
517   
518   if (!fRunLoader) {
519       AliError("No run loader is found !"); 
520       return kFALSE;
521     }
522     // read run number from gAlice
523     if(fRunLoader->GetAliRun())
524       AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
525     else {
526       if(fRawReader) {
527         if(fRawReader->NextEvent()) {
528           AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
529           fRawReader->RewindEvents();
530         }
531         else {
532           if(man->GetRun() > 0) {
533             AliWarning("No raw events is found ! Using settings in AliCDBManager !");
534             man->Print();  
535             return kTRUE;
536           }
537           else {
538             AliWarning("Neither raw events nor settings in AliCDBManager are found !");
539             return kFALSE;
540           }
541         }
542       }
543       else {
544         AliError("Neither gAlice nor RawReader objects are found !");
545         return kFALSE;
546       }
547   }
548
549   man->Print();  
550   
551   return kTRUE;
552 }
553
554 //_____________________________________________________________________________
555 void AliReconstruction::SetCDBLock() {
556   // Set CDB lock: from now on it is forbidden to reset the run number
557   // or the default storage or to activate any further storage!
558   
559   AliCDBManager::Instance()->SetLock(1);
560 }
561
562 //_____________________________________________________________________________
563 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
564 {
565   // Read the alignment objects from CDB.
566   // Each detector is supposed to have the
567   // alignment objects in DET/Align/Data CDB path.
568   // All the detector objects are then collected,
569   // sorted by geometry level (starting from ALIC) and
570   // then applied to the TGeo geometry.
571   // Finally an overlaps check is performed.
572
573   // Load alignment data from CDB and fill fAlignObjArray 
574   if(fLoadAlignFromCDB){
575         
576     TString detStr = detectors;
577     TString loadAlObjsListOfDets = "";
578     
579     for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
580       if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
581       loadAlObjsListOfDets += fgkDetectorName[iDet];
582       loadAlObjsListOfDets += " ";
583     } // end loop over detectors
584     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
585     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
586   }else{
587     // Check if the array with alignment objects was
588     // provided by the user. If yes, apply the objects
589     // to the present TGeo geometry
590     if (fAlignObjArray) {
591       if (gGeoManager && gGeoManager->IsClosed()) {
592         if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
593           AliError("The misalignment of one or more volumes failed!"
594                    "Compare the list of simulated detectors and the list of detector alignment data!");
595           return kFALSE;
596         }
597       }
598       else {
599         AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
600         return kFALSE;
601       }
602     }
603   }
604   
605   delete fAlignObjArray; fAlignObjArray=0;
606
607   return kTRUE;
608 }
609
610 //_____________________________________________________________________________
611 void AliReconstruction::SetGAliceFile(const char* fileName)
612 {
613 // set the name of the galice file
614
615   fGAliceFileName = fileName;
616 }
617
618 //_____________________________________________________________________________
619 void AliReconstruction::SetInput(const char* input) 
620 {
621   // In case the input string starts with 'mem://', we run in an online mode
622   // and AliRawReaderDateOnline object is created. In all other cases a raw-data
623   // file is assumed. One can give as an input:
624   // mem://: - events taken from DAQ monitoring libs online
625   //  or
626   // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
627   fInput = input;
628 }
629
630 //_____________________________________________________________________________
631 void AliReconstruction::SetOption(const char* detector, const char* option)
632 {
633 // set options for the reconstruction of a detector
634
635   TObject* obj = fOptions.FindObject(detector);
636   if (obj) fOptions.Remove(obj);
637   fOptions.Add(new TNamed(detector, option));
638 }
639
640 //_____________________________________________________________________________
641 Bool_t AliReconstruction::Run(const char* input)
642 {
643   // Run Run Run
644   AliCodeTimerAuto("");
645
646   if (!InitRun(input)) return kFALSE;
647
648   //******* The loop over events
649   Int_t iEvent = 0;
650   while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
651          (fRawReader && fRawReader->NextEvent())) {
652     if (!RunEvent(iEvent)) return kFALSE;
653     iEvent++;
654   }
655
656   if (!FinishRun()) return kFALSE;
657
658   return kTRUE;
659 }
660
661 //_____________________________________________________________________________
662 Bool_t AliReconstruction::InitRun(const char* input)
663 {
664   // Initialize all the stuff before
665   // going into the event loop
666   // If the second argument is given, the first one is ignored and
667   // the reconstruction works in an online mode
668   AliCodeTimerAuto("");
669
670   // Overwrite the previous setting
671   if (input) fInput = input;
672
673   // set the input in case of raw data
674   fRawReader = AliRawReader::Create(fInput.Data());
675   if (!fRawReader)
676     AliInfo("Reconstruction will run over digits");
677
678   if (!fEquipIdMap.IsNull() && fRawReader)
679     fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
680
681   if (!fUseHLTData.IsNull()) {
682     // create the RawReaderHLT which performs redirection of HLT input data for
683     // the specified detectors
684     AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
685     if (pRawReader) {
686       fParentRawReader=fRawReader;
687       fRawReader=pRawReader;
688     } else {
689       AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
690     }
691   }
692
693    AliSysInfo::AddStamp("Start");
694   // get the run loader
695   if (!InitRunLoader()) return kFALSE;
696    AliSysInfo::AddStamp("LoadLoader");
697
698   // Initialize the CDB storage
699   InitCDB();
700   
701   AliSysInfo::AddStamp("LoadCDB");
702
703   // Set run number in CDBManager (if it is not already set by the user)
704   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
705   
706   // Set CDB lock: from now on it is forbidden to reset the run number
707   // or the default storage or to activate any further storage!
708   SetCDBLock();
709   
710   // Import ideal TGeo geometry and apply misalignment
711   if (!gGeoManager) {
712     TString geom(gSystem->DirName(fGAliceFileName));
713     geom += "/geometry.root";
714     AliGeomManager::LoadGeometry(geom.Data());
715     if (!gGeoManager) if (fStopOnError) return kFALSE;
716   }
717
718   if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
719    AliSysInfo::AddStamp("LoadGeom");
720
721   //QA
722   AliQADataMakerSteer qas ; 
723   if (fRunQA && fRawReader) { 
724     qas.Run(fRunLocalReconstruction, fRawReader) ; 
725         fSameQACycle = kTRUE ; 
726   }
727   // checking the QA of previous steps
728   //CheckQA() ; 
729  
730   /*
731   // local reconstruction
732   if (!fRunLocalReconstruction.IsNull()) {
733     if (!RunLocalReconstruction(fRunLocalReconstruction)) {
734       if (fStopOnError) {CleanUp(); return kFALSE;}
735     }
736   }
737   */
738
739   // get vertexer
740   if (fRunVertexFinder && !CreateVertexer()) {
741     if (fStopOnError) {
742       CleanUp(); 
743       return kFALSE;
744     }
745   }
746    AliSysInfo::AddStamp("Vertexer");
747
748   // get trackers
749   if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
750     if (fStopOnError) {
751       CleanUp(); 
752       return kFALSE;
753     }      
754   }
755    AliSysInfo::AddStamp("LoadTrackers");
756
757   // get the possibly already existing ESD file and tree
758   fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
759   if (!gSystem->AccessPathName("AliESDs.root")){
760     gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
761     ffileOld = TFile::Open("AliESDs.old.root");
762     if (ffileOld && ffileOld->IsOpen()) {
763       ftreeOld = (TTree*) ffileOld->Get("esdTree");
764       if (ftreeOld)fesd->ReadFromTree(ftreeOld);
765       fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
766       if (fhlttreeOld)  fhltesd->ReadFromTree(fhlttreeOld);
767     }
768   }
769
770   // create the ESD output file and tree
771   ffile = TFile::Open("AliESDs.root", "RECREATE");
772   ffile->SetCompressionLevel(2);
773   if (!ffile->IsOpen()) {
774     AliError("opening AliESDs.root failed");
775     if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}    
776   }
777
778   ftree = new TTree("esdTree", "Tree with ESD objects");
779   fesd = new AliESDEvent();
780   fesd->CreateStdContent();
781   fesd->WriteToTree(ftree);
782
783   fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
784   fhltesd = new AliESDEvent();
785   fhltesd->CreateStdContent();
786   fhltesd->WriteToTree(fhlttree);
787
788   /* CKB Why?
789   delete esd; delete hltesd;
790   esd = NULL; hltesd = NULL;
791   */
792   // create the branch with ESD additions
793
794
795
796   if (fWriteESDfriend) {
797     fesdf = new AliESDfriend();
798     TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
799     br->SetFile("AliESDfriends.root");
800     fesd->AddObject(fesdf);
801   }
802
803   // Get the GRP CDB entry
804   AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
805         
806   if (entryGRP) 
807         fGRPData = dynamic_cast<TMap*> (entryGRP->GetObject());  
808   
809   if (!fGRPData)
810         AliError("No GRP entry found in OCDB!");
811
812   // Get the diamond profile from OCDB
813   AliCDBEntry* entry = AliCDBManager::Instance()
814         ->Get("GRP/Calib/MeanVertex");
815         
816   if(entry) {
817         fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());  
818   } else {
819         AliError("No diamond profile found in OCDB!");
820   }
821
822   entry = 0;
823   entry = AliCDBManager::Instance()
824         ->Get("GRP/Calib/MeanVertexTPC");
825         
826   if(entry) {
827         fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());  
828   } else {
829         AliError("No diamond profile found in OCDB!");
830   }
831
832   ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
833   if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
834
835   if (fRawReader) fRawReader->RewindEvents();
836
837   ProcInfo_t ProcInfo;
838   gSystem->GetProcInfo(&ProcInfo);
839   AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
840   
841
842   //Initialize the QA and start of cycle for out-of-cycle QA
843   if (fRunQA) {
844      TString detStr(fFillESD); 
845      for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
846         if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
847         AliQADataMakerRec *qadm = GetQADataMaker(iDet);  
848         if (!qadm) continue;
849         AliInfo(Form("Initializing the QA data maker for %s", 
850                fgkDetectorName[iDet]));
851         qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
852         qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
853         if (!fInLoopQA) {
854                         qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
855                         qadm->StartOfCycle(AliQA::kESDS,"same");
856         }
857      }
858           if (fRunGlobalQA) {
859                   AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
860                   AliInfo(Form("Initializing the global QA data maker"));
861                   TObjArray *arr=
862                         qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
863                   AliTracker::SetResidualsArray(arr);
864                   qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
865                   if (!fInLoopQA) {
866                           qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
867                           qadm->StartOfCycle(AliQA::kESDS, "same");
868                   }
869           }
870           if (!fInLoopQA) 
871                   fSameQACycle = kTRUE; 
872   }
873
874   //Initialize the Plane Efficiency framework
875   if (fRunPlaneEff && !InitPlaneEff()) {
876     if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
877   }
878
879   if (strcmp(gProgName,"alieve") == 0)
880     fRunAliEVE = InitAliEVE();
881
882   return kTRUE;
883 }
884
885 //_____________________________________________________________________________
886 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
887 {
888   // run the reconstruction over a single event
889   // The event loop is steered in Run method
890
891   AliCodeTimerAuto("");
892
893   if (iEvent >= fRunLoader->GetNumberOfEvents()) {
894     fRunLoader->SetEventNumber(iEvent);
895     fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(), 
896                                    iEvent, iEvent);
897     //??      fRunLoader->MakeTree("H");
898     fRunLoader->TreeE()->Fill();
899   }
900
901   if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
902     // copy old ESD to the new one
903     if (ftreeOld) {
904       fesd->ReadFromTree(ftreeOld);
905       ftreeOld->GetEntry(iEvent);
906       ftree->Fill();
907     }
908     if (fhlttreeOld) {
909       fhltesd->ReadFromTree(fhlttreeOld);
910       fhlttreeOld->GetEntry(iEvent);
911       fhlttree->Fill();
912     }
913     return kTRUE;
914   }
915
916   AliInfo(Form("processing event %d", iEvent));
917
918     //Start of cycle for the in-loop QA
919     if (fInLoopQA) {
920        if (fRunQA) {
921           TString detStr(fFillESD); 
922           for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
923              if (!IsSelected(fgkDetectorName[iDet], detStr)) 
924                                  continue;
925              AliQADataMakerRec *qadm = GetQADataMaker(iDet);  
926              if (!qadm) 
927                                  continue;
928              qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
929              qadm->StartOfCycle(AliQA::kESDS, "same") ;         
930           }
931                    if (fRunGlobalQA) {
932                            AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
933                            qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
934                            qadm->StartOfCycle(AliQA::kESDS, "same");
935                    }               
936            }
937     }
938
939     fRunLoader->GetEvent(iEvent);
940
941     char aFileName[256];
942     sprintf(aFileName, "ESD_%d.%d_final.root", 
943             fRunLoader->GetHeader()->GetRun(), 
944             fRunLoader->GetHeader()->GetEventNrInRun());
945     if (!gSystem->AccessPathName(aFileName)) return kTRUE;
946
947     // local single event reconstruction
948     if (!fRunLocalReconstruction.IsNull()) {
949       TString detectors=fRunLocalReconstruction;
950       // run HLT event reconstruction first
951       // ;-( IsSelected changes the string
952       if (IsSelected("HLT", detectors) &&
953           !RunLocalEventReconstruction("HLT")) {
954         if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
955       }
956       detectors=fRunLocalReconstruction;
957       detectors.ReplaceAll("HLT", "");
958       if (!RunLocalEventReconstruction(detectors)) {
959         if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
960       }
961     }
962
963     fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
964     fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
965     fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
966     fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
967     
968     // Set magnetic field from the tracker
969     fesd->SetMagneticField(AliTracker::GetBz());
970     fhltesd->SetMagneticField(AliTracker::GetBz());
971
972     
973     
974     // Fill raw-data error log into the ESD
975     if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
976
977     // vertex finder
978     if (fRunVertexFinder) {
979       if (!ReadESD(fesd, "vertex")) {
980         if (!RunVertexFinder(fesd)) {
981           if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
982         }
983         if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
984       }
985     }
986
987     // Muon tracking
988     if (!fRunTracking.IsNull()) {
989       if (fRunMuonTracking) {
990         if (!RunMuonTracking(fesd)) {
991           if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
992         }
993       }
994     }
995
996     // barrel tracking
997     if (!fRunTracking.IsNull()) {
998       if (!ReadESD(fesd, "tracking")) {
999         if (!RunTracking(fesd)) {
1000           if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1001         }
1002         if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1003       }
1004     }
1005
1006     // fill ESD
1007     if (!fFillESD.IsNull()) {
1008       TString detectors=fFillESD;
1009       // run HLT first and on hltesd
1010       // ;-( IsSelected changes the string
1011       if (IsSelected("HLT", detectors) &&
1012           !FillESD(fhltesd, "HLT")) {
1013         if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1014       }
1015       detectors=fFillESD;
1016       detectors.ReplaceAll("HLT", "");
1017       if (!FillESD(fesd, detectors)) {
1018         if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1019       }
1020     }
1021   
1022     // fill Event header information from the RawEventHeader
1023     if (fRawReader){FillRawEventHeaderESD(fesd);}
1024
1025     // combined PID
1026     AliESDpid::MakePID(fesd);
1027     if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1028
1029     if (fFillTriggerESD) {
1030       if (!ReadESD(fesd, "trigger")) {
1031         if (!FillTriggerESD(fesd)) {
1032           if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1033         }
1034         if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1035       }
1036     }
1037
1038     ffile->cd();
1039
1040     //
1041     // Propagate track to the beam pipe  (if not already done by ITS)
1042     //
1043     const Int_t ntracks = fesd->GetNumberOfTracks();
1044     const Double_t kBz = fesd->GetMagneticField();
1045     const Double_t kRadius  = 2.8; //something less than the beam pipe radius
1046
1047     TObjArray trkArray;
1048     UShort_t *selectedIdx=new UShort_t[ntracks];
1049
1050     for (Int_t itrack=0; itrack<ntracks; itrack++){
1051       const Double_t kMaxStep = 5;   //max step over the material
1052       Bool_t ok;
1053
1054       AliESDtrack *track = fesd->GetTrack(itrack);
1055       if (!track) continue;
1056
1057       AliExternalTrackParam *tpcTrack =
1058            (AliExternalTrackParam *)track->GetTPCInnerParam();
1059       ok = kFALSE;
1060       if (tpcTrack)
1061         ok = AliTracker::
1062           PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1063
1064
1065
1066       if (ok) {
1067         Int_t n=trkArray.GetEntriesFast();
1068         selectedIdx[n]=track->GetID();
1069         trkArray.AddLast(tpcTrack);
1070       }
1071
1072       if (track->GetX() < kRadius) continue;
1073
1074       ok = AliTracker::
1075            PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1076       if (ok) {
1077          track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kRadius);
1078       }
1079     }
1080
1081     //
1082     // Improve the reconstructed primary vertex position using the tracks
1083     //
1084     TObject *obj = fOptions.FindObject("ITS");
1085     if (obj) {
1086       TString optITS = obj->GetTitle();
1087       if (optITS.Contains("cosmics") || optITS.Contains("COSMICS")) 
1088         fRunVertexFinderTracks=kFALSE;
1089     }
1090     if (fRunVertexFinderTracks) {
1091        // TPC + ITS primary vertex
1092        ftVertexer->SetITSrefitRequired();
1093        if(fDiamondProfile && fMeanVertexConstraint) {
1094          ftVertexer->SetVtxStart(fDiamondProfile);
1095        } else {
1096          ftVertexer->SetConstraintOff();
1097        }
1098        AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1099        if (pvtx) {
1100           if (pvtx->GetStatus()) {
1101              fesd->SetPrimaryVertex(pvtx);
1102              for (Int_t i=0; i<ntracks; i++) {
1103                  AliESDtrack *t = fesd->GetTrack(i);
1104                  t->RelateToVertex(pvtx, kBz, kRadius);
1105              } 
1106           }
1107        }
1108
1109        // TPC-only primary vertex
1110        ftVertexer->SetITSrefitNotRequired();
1111        if(fDiamondProfileTPC && fMeanVertexConstraint) {
1112          ftVertexer->SetVtxStart(fDiamondProfileTPC);
1113        } else {
1114          ftVertexer->SetConstraintOff();
1115        }
1116        pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1117        if (pvtx) {
1118           if (pvtx->GetStatus()) {
1119              fesd->SetPrimaryVertexTPC(pvtx);
1120              Int_t nsel=trkArray.GetEntriesFast();
1121              for (Int_t i=0; i<nsel; i++) {
1122                  AliExternalTrackParam *t = 
1123                    (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1124                  t->PropagateToDCA(pvtx, kBz, kRadius);
1125              } 
1126           }
1127        }
1128
1129     }
1130     delete[] selectedIdx;
1131
1132     if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1133     
1134
1135     if (fRunV0Finder) {
1136        // V0 finding
1137        AliV0vertexer vtxer;
1138        vtxer.Tracks2V0vertices(fesd);
1139
1140        if (fRunCascadeFinder) {
1141           // Cascade finding
1142           AliCascadeVertexer cvtxer;
1143           cvtxer.V0sTracks2CascadeVertices(fesd);
1144        }
1145     }
1146  
1147     // write ESD
1148     if (fCleanESD) CleanESD(fesd);
1149
1150         if (fRunQA) {
1151                 if (fRunGlobalQA) {
1152                         AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1153                         if (qadm) qadm->Exec(AliQA::kESDS, fesd);
1154                 }
1155         }
1156
1157     if (fWriteESDfriend) {
1158       fesdf->~AliESDfriend();
1159       new (fesdf) AliESDfriend(); // Reset...
1160       fesd->GetESDfriend(fesdf);
1161     }
1162     ftree->Fill();
1163
1164     // write HLT ESD
1165     fhlttree->Fill();
1166
1167     // call AliEVE
1168     if (fRunAliEVE) RunAliEVE();
1169
1170     if (fCheckPointLevel > 0)  WriteESD(fesd, "final"); 
1171     fesd->Reset();
1172     fhltesd->Reset();
1173     if (fWriteESDfriend) {
1174       fesdf->~AliESDfriend();
1175       new (fesdf) AliESDfriend(); // Reset...
1176     }
1177  
1178     ProcInfo_t ProcInfo;
1179     gSystem->GetProcInfo(&ProcInfo);
1180     AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1181   
1182
1183   // End of cycle for the in-loop QA
1184      if (fInLoopQA) {
1185         if (fRunQA) {
1186            RunQA(fFillESD.Data(), fesd);
1187            TString detStr(fFillESD); 
1188            for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1189                            if (!IsSelected(fgkDetectorName[iDet], detStr)) 
1190                                    continue;
1191                            AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1192                            if (!qadm) 
1193                                    continue;
1194                            qadm->EndOfCycle(AliQA::kRECPOINTS);
1195                            qadm->EndOfCycle(AliQA::kESDS);
1196                            qadm->Finish();
1197                    }
1198         }
1199         if (fRunGlobalQA) {
1200            AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1201            if (qadm) {
1202               qadm->EndOfCycle(AliQA::kRECPOINTS);
1203               qadm->EndOfCycle(AliQA::kESDS);
1204               qadm->Finish();
1205            }
1206         }
1207      }
1208
1209      return kTRUE;
1210 }
1211
1212 //_____________________________________________________________________________
1213 Bool_t AliReconstruction::FinishRun()
1214 {
1215   // Finalize the run
1216   // Called after the exit
1217   // from the event loop
1218   AliCodeTimerAuto("");
1219
1220   if (fIsNewRunLoader) { // galice.root didn't exist
1221     fRunLoader->WriteHeader("OVERWRITE");
1222     fRunLoader->CdGAFile();
1223     fRunLoader->Write(0, TObject::kOverwrite);
1224   }
1225
1226   ftree->GetUserInfo()->Add(fesd);
1227   fhlttree->GetUserInfo()->Add(fhltesd);
1228   
1229   const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();       
1230   const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();   
1231                  
1232    TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());    
1233    cdbMapCopy->SetOwner(1);      
1234    cdbMapCopy->SetName("cdbMap");        
1235    TIter iter(cdbMap->GetTable());       
1236          
1237    TPair* pair = 0;      
1238    while((pair = dynamic_cast<TPair*> (iter.Next()))){   
1239          TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());   
1240          TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());         
1241          cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));  
1242    }     
1243          
1244    TList *cdbListCopy = new TList();     
1245    cdbListCopy->SetOwner(1);     
1246    cdbListCopy->SetName("cdbList");      
1247          
1248    TIter iter2(cdbList);         
1249          
1250    AliCDBId* id=0;       
1251    while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){         
1252          cdbListCopy->Add(new TObjString(id->ToString().Data()));        
1253    }     
1254          
1255    ftree->GetUserInfo()->Add(cdbMapCopy);        
1256    ftree->GetUserInfo()->Add(cdbListCopy);
1257
1258
1259   if(fESDPar.Contains("ESD.par")){
1260     AliInfo("Attaching ESD.par to Tree");
1261     TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1262     ftree->GetUserInfo()->Add(fn);
1263   }
1264
1265
1266   ffile->cd();
1267
1268   if (fWriteESDfriend)
1269     ftree->SetBranchStatus("ESDfriend*",0);
1270   // we want to have only one tree version number
1271   ftree->Write(ftree->GetName(),TObject::kOverwrite);
1272   fhlttree->Write();
1273
1274 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1275   if (fRunPlaneEff && !FinishPlaneEff()) {
1276    AliWarning("Finish PlaneEff evaluation failed");
1277   }
1278
1279   gROOT->cd();
1280   CleanUp(ffile, ffileOld);
1281     
1282   if (fWriteAOD) {
1283     AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1284   }
1285
1286   // Create tags for the events in the ESD tree (the ESD tree is always present)
1287   // In case of empty events the tags will contain dummy values
1288   AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1289   esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1290   if (fWriteAOD) {
1291     AliWarning("AOD tag creation not supported anymore during reconstruction.");
1292   }
1293
1294   //Finish QA and end of cycle for out-of-loop QA
1295   if (!fInLoopQA) {
1296           if (fRunQA) {
1297                   AliQADataMakerSteer qas;
1298                   qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1299                   //qas.Reset() ;
1300                   qas.Run(fRunTracking.Data(), AliQA::kESDS, fSameQACycle);
1301                   if (fRunGlobalQA) {
1302                          AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1303                           if (qadm) {
1304                                   qadm->EndOfCycle(AliQA::kRECPOINTS);
1305                                   qadm->EndOfCycle(AliQA::kESDS);
1306                                   qadm->Finish();
1307                           }
1308                   }
1309           }
1310   }
1311   
1312   // Cleanup of CDB manager: cache and active storages!
1313   AliCDBManager::Instance()->ClearCache();
1314   
1315   return kTRUE;
1316 }
1317
1318
1319 //_____________________________________________________________________________
1320 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1321 {
1322 // run the local reconstruction
1323   static Int_t eventNr=0;
1324   AliCodeTimerAuto("")
1325
1326  //  AliCDBManager* man = AliCDBManager::Instance();
1327 //   Bool_t origCache = man->GetCacheFlag();
1328
1329 //   TString detStr = detectors;
1330 //   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1331 //     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1332 //     AliReconstructor* reconstructor = GetReconstructor(iDet);
1333 //     if (!reconstructor) continue;
1334 //     if (reconstructor->HasLocalReconstruction()) continue;
1335
1336 //     AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1337 //     AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1338     
1339 //     AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));                          
1340 //     AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1341
1342 //     man->SetCacheFlag(kTRUE);
1343 //     TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1344 //     man->GetAll(calibPath); // entries are cached!
1345
1346 //     AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1347      
1348 //     if (fRawReader) {
1349 //       fRawReader->RewindEvents();
1350 //       reconstructor->Reconstruct(fRunLoader, fRawReader);
1351 //     } else {
1352 //       reconstructor->Reconstruct(fRunLoader);
1353 //     }
1354      
1355 //      AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1356     // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1357
1358 //     // unload calibration data
1359 //     man->UnloadFromCache(calibPath);
1360 //     //man->ClearCache();
1361 //   }
1362
1363 //   man->SetCacheFlag(origCache);
1364
1365 //   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1366 //     AliError(Form("the following detectors were not found: %s",
1367 //                   detStr.Data()));
1368 //     if (fStopOnError) return kFALSE;
1369 //   }
1370
1371           eventNr++;
1372   return kTRUE;
1373 }
1374
1375 //_____________________________________________________________________________
1376 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1377 {
1378 // run the local reconstruction
1379
1380   static Int_t eventNr=0;
1381   AliCodeTimerAuto("")
1382
1383   TString detStr = detectors;
1384   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1385     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1386     AliReconstructor* reconstructor = GetReconstructor(iDet);
1387     if (!reconstructor) continue;
1388     AliLoader* loader = fLoader[iDet];
1389     // Matthias April 2008: temporary fix to run HLT reconstruction
1390     // although the HLT loader is missing
1391     if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1392       if (fRawReader) {
1393         reconstructor->Reconstruct(fRawReader, NULL);
1394       } else {
1395         TTree* dummy=NULL;
1396         reconstructor->Reconstruct(dummy, NULL);
1397       }
1398       continue;
1399     }
1400     if (!loader) {
1401       AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1402       continue;
1403     }
1404     // conversion of digits
1405     if (fRawReader && reconstructor->HasDigitConversion()) {
1406       AliInfo(Form("converting raw data digits into root objects for %s", 
1407                    fgkDetectorName[iDet]));
1408       AliCodeTimerAuto(Form("converting raw data digits into root objects for %s", 
1409                             fgkDetectorName[iDet]));
1410       loader->LoadDigits("update");
1411       loader->CleanDigits();
1412       loader->MakeDigitsContainer();
1413       TTree* digitsTree = loader->TreeD();
1414       reconstructor->ConvertDigits(fRawReader, digitsTree);
1415       loader->WriteDigits("OVERWRITE");
1416       loader->UnloadDigits();
1417     }
1418     // local reconstruction
1419     AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1420     AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1421     loader->LoadRecPoints("update");
1422     loader->CleanRecPoints();
1423     loader->MakeRecPointsContainer();
1424     TTree* clustersTree = loader->TreeR();
1425     if (fRawReader && !reconstructor->HasDigitConversion()) {
1426       reconstructor->Reconstruct(fRawReader, clustersTree);
1427     } else {
1428       loader->LoadDigits("read");
1429       TTree* digitsTree = loader->TreeD();
1430       if (!digitsTree) {
1431         AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1432         if (fStopOnError) return kFALSE;
1433       } else {
1434         reconstructor->Reconstruct(digitsTree, clustersTree);
1435       }
1436       loader->UnloadDigits();
1437     }
1438
1439     // In-loop QA for local reconstrucion 
1440     if (fRunQA && fInLoopQA) {
1441        AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1442        if (qadm) {
1443           //AliCodeTimerStart
1444           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1445           //AliInfo
1446           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1447
1448           qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1449  
1450           //AliCodeTimerStop
1451           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1452        }
1453     }
1454
1455     loader->WriteRecPoints("OVERWRITE");
1456     loader->UnloadRecPoints();
1457     AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1458   }
1459
1460   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1461     AliError(Form("the following detectors were not found: %s",
1462                   detStr.Data()));
1463     if (fStopOnError) return kFALSE;
1464   }
1465   eventNr++;
1466   return kTRUE;
1467 }
1468
1469 //_____________________________________________________________________________
1470 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1471 {
1472 // run the barrel tracking
1473
1474   AliCodeTimerAuto("")
1475
1476   AliESDVertex* vertex = NULL;
1477   Double_t vtxPos[3] = {0, 0, 0};
1478   Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1479   TArrayF mcVertex(3); 
1480   if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1481     fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1482     for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1483   }
1484
1485   if (fVertexer) {
1486     if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1487     AliInfo("running the ITS vertex finder");
1488     if (fLoader[0]) fLoader[0]->LoadRecPoints();
1489     vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1490     if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1491     if(!vertex){
1492       AliWarning("Vertex not found");
1493       vertex = new AliESDVertex();
1494       vertex->SetName("default");
1495     }
1496     else {
1497       vertex->SetName("reconstructed");
1498     }
1499
1500   } else {
1501     AliInfo("getting the primary vertex from MC");
1502     vertex = new AliESDVertex(vtxPos, vtxErr);
1503   }
1504
1505   if (vertex) {
1506     vertex->GetXYZ(vtxPos);
1507     vertex->GetSigmaXYZ(vtxErr);
1508   } else {
1509     AliWarning("no vertex reconstructed");
1510     vertex = new AliESDVertex(vtxPos, vtxErr);
1511   }
1512   esd->SetPrimaryVertexSPD(vertex);
1513   // if SPD multiplicity has been determined, it is stored in the ESD
1514   AliMultiplicity *mult = fVertexer->GetMultiplicity();
1515   if(mult)esd->SetMultiplicity(mult);
1516
1517   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1518     if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1519   }  
1520   delete vertex;
1521
1522   return kTRUE;
1523 }
1524
1525 //_____________________________________________________________________________
1526 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1527 {
1528 // run the HLT barrel tracking
1529
1530   AliCodeTimerAuto("")
1531
1532   if (!fRunLoader) {
1533     AliError("Missing runLoader!");
1534     return kFALSE;
1535   }
1536
1537   AliInfo("running HLT tracking");
1538
1539   // Get a pointer to the HLT reconstructor
1540   AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1541   if (!reconstructor) return kFALSE;
1542
1543   // TPC + ITS
1544   for (Int_t iDet = 1; iDet >= 0; iDet--) {
1545     TString detName = fgkDetectorName[iDet];
1546     AliDebug(1, Form("%s HLT tracking", detName.Data()));
1547     reconstructor->SetOption(detName.Data());
1548     AliTracker *tracker = reconstructor->CreateTracker();
1549     if (!tracker) {
1550       AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1551       if (fStopOnError) return kFALSE;
1552       continue;
1553     }
1554     Double_t vtxPos[3];
1555     Double_t vtxErr[3]={0.005,0.005,0.010};
1556     const AliESDVertex *vertex = esd->GetVertex();
1557     vertex->GetXYZ(vtxPos);
1558     tracker->SetVertex(vtxPos,vtxErr);
1559     if(iDet != 1) {
1560       fLoader[iDet]->LoadRecPoints("read");
1561       TTree* tree = fLoader[iDet]->TreeR();
1562       if (!tree) {
1563         AliError(Form("Can't get the %s cluster tree", detName.Data()));
1564         return kFALSE;
1565       }
1566       tracker->LoadClusters(tree);
1567     }
1568     if (tracker->Clusters2Tracks(esd) != 0) {
1569       AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1570       return kFALSE;
1571     }
1572     if(iDet != 1) {
1573       tracker->UnloadClusters();
1574     }
1575     delete tracker;
1576   }
1577
1578   return kTRUE;
1579 }
1580
1581 //_____________________________________________________________________________
1582 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1583 {
1584 // run the muon spectrometer tracking
1585
1586   AliCodeTimerAuto("")
1587
1588   if (!fRunLoader) {
1589     AliError("Missing runLoader!");
1590     return kFALSE;
1591   }
1592   Int_t iDet = 7; // for MUON
1593
1594   AliInfo("is running...");
1595
1596   // Get a pointer to the MUON reconstructor
1597   AliReconstructor *reconstructor = GetReconstructor(iDet);
1598   if (!reconstructor) return kFALSE;
1599
1600   
1601   TString detName = fgkDetectorName[iDet];
1602   AliDebug(1, Form("%s tracking", detName.Data()));
1603   AliTracker *tracker =  reconstructor->CreateTracker();
1604   if (!tracker) {
1605     AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1606     return kFALSE;
1607   }
1608      
1609   // read RecPoints
1610   fLoader[iDet]->LoadRecPoints("read");  
1611
1612   tracker->LoadClusters(fLoader[iDet]->TreeR());
1613   
1614   Int_t rv = tracker->Clusters2Tracks(esd);
1615   
1616   if ( rv )
1617   {
1618     AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1619     return kFALSE;
1620   }
1621   
1622   fLoader[iDet]->UnloadRecPoints();
1623
1624   tracker->UnloadClusters();
1625   
1626   delete tracker;
1627   
1628   return kTRUE;
1629 }
1630
1631
1632 //_____________________________________________________________________________
1633 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1634 {
1635 // run the barrel tracking
1636   static Int_t eventNr=0;
1637   AliCodeTimerAuto("")
1638
1639   AliInfo("running tracking");
1640
1641   //Fill the ESD with the T0 info (will be used by the TOF) 
1642   if (fReconstructor[11] && fLoader[11]) {
1643     fLoader[11]->LoadRecPoints("READ");
1644     TTree *treeR = fLoader[11]->TreeR();
1645     GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1646   }
1647
1648   // pass 1: TPC + ITS inwards
1649   for (Int_t iDet = 1; iDet >= 0; iDet--) {
1650     if (!fTracker[iDet]) continue;
1651     AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1652
1653     // load clusters
1654     fLoader[iDet]->LoadRecPoints("read");
1655     AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1656     TTree* tree = fLoader[iDet]->TreeR();
1657     if (!tree) {
1658       AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1659       return kFALSE;
1660     }
1661     fTracker[iDet]->LoadClusters(tree);
1662     AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1663     // run tracking
1664     if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1665       AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1666       return kFALSE;
1667     }
1668     if (fCheckPointLevel > 1) {
1669       WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1670     }
1671     // preliminary PID in TPC needed by the ITS tracker
1672     if (iDet == 1) {
1673       GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1674       AliESDpid::MakePID(esd);
1675     } 
1676     AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1677   }
1678
1679   // pass 2: ALL backwards
1680
1681   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1682     if (!fTracker[iDet]) continue;
1683     AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1684
1685     // load clusters
1686     if (iDet > 1) {     // all except ITS, TPC
1687       TTree* tree = NULL;
1688       fLoader[iDet]->LoadRecPoints("read");
1689       AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1690       tree = fLoader[iDet]->TreeR();
1691       if (!tree) {
1692         AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1693         return kFALSE;
1694       }
1695       fTracker[iDet]->LoadClusters(tree); 
1696       AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1697     }
1698
1699     // run tracking
1700     if (iDet>1) // start filling residuals for the "outer" detectors
1701     if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);     
1702
1703     if (fTracker[iDet]->PropagateBack(esd) != 0) {
1704       AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1705       //      return kFALSE;
1706     }
1707     if (fCheckPointLevel > 1) {
1708       WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1709     }
1710
1711     // unload clusters
1712     if (iDet > 2) {     // all except ITS, TPC, TRD
1713       fTracker[iDet]->UnloadClusters();
1714       fLoader[iDet]->UnloadRecPoints();
1715     }
1716     // updated PID in TPC needed by the ITS tracker -MI
1717     if (iDet == 1) {
1718       GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1719       AliESDpid::MakePID(esd);
1720     }
1721     AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1722   }
1723   //stop filling residuals for the "outer" detectors
1724   if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);     
1725
1726   // write space-points to the ESD in case alignment data output
1727   // is switched on
1728   if (fWriteAlignmentData)
1729     WriteAlignmentData(esd);
1730
1731   // pass 3: TRD + TPC + ITS refit inwards
1732
1733   for (Int_t iDet = 2; iDet >= 0; iDet--) {
1734     if (!fTracker[iDet]) continue;
1735     AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1736
1737     // run tracking
1738     if (iDet<2) // start filling residuals for TPC and ITS
1739     if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);     
1740
1741     if (fTracker[iDet]->RefitInward(esd) != 0) {
1742       AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1743       //      return kFALSE;
1744     }
1745     // run postprocessing
1746     if (fTracker[iDet]->PostProcess(esd) != 0) {
1747       AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1748       //      return kFALSE;
1749     }
1750     if (fCheckPointLevel > 1) {
1751       WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1752     }
1753     AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1754     // unload clusters
1755     fTracker[iDet]->UnloadClusters();
1756     AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1757     fLoader[iDet]->UnloadRecPoints();
1758     AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1759   }
1760   // stop filling residuals for TPC and ITS
1761   if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);     
1762
1763   eventNr++;
1764   return kTRUE;
1765 }
1766
1767 //_____________________________________________________________________________
1768 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1769   //
1770   // Remove the data which are not needed for the physics analysis.
1771   //
1772
1773   Int_t nTracks=esd->GetNumberOfTracks();
1774   Int_t nV0s=esd->GetNumberOfV0s();
1775   AliInfo
1776   (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1777
1778   Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1779   Bool_t rc=esd->Clean(cleanPars);
1780
1781   nTracks=esd->GetNumberOfTracks();
1782   nV0s=esd->GetNumberOfV0s();
1783   AliInfo
1784   (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1785
1786   return rc;
1787 }
1788
1789 //_____________________________________________________________________________
1790 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1791 {
1792 // fill the event summary data
1793
1794   AliCodeTimerAuto("")
1795     static Int_t eventNr=0; 
1796   TString detStr = detectors;
1797   
1798   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1799   if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1800     AliReconstructor* reconstructor = GetReconstructor(iDet);
1801     if (!reconstructor) continue;
1802     if (!ReadESD(esd, fgkDetectorName[iDet])) {
1803       AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1804       TTree* clustersTree = NULL;
1805       if (fLoader[iDet]) {
1806         fLoader[iDet]->LoadRecPoints("read");
1807         clustersTree = fLoader[iDet]->TreeR();
1808         if (!clustersTree) {
1809           AliError(Form("Can't get the %s clusters tree", 
1810                         fgkDetectorName[iDet]));
1811           if (fStopOnError) return kFALSE;
1812         }
1813       }
1814       if (fRawReader && !reconstructor->HasDigitConversion()) {
1815         reconstructor->FillESD(fRawReader, clustersTree, esd);
1816       } else {
1817         TTree* digitsTree = NULL;
1818         if (fLoader[iDet]) {
1819           fLoader[iDet]->LoadDigits("read");
1820           digitsTree = fLoader[iDet]->TreeD();
1821           if (!digitsTree) {
1822             AliError(Form("Can't get the %s digits tree", 
1823                           fgkDetectorName[iDet]));
1824             if (fStopOnError) return kFALSE;
1825           }
1826         }
1827         reconstructor->FillESD(digitsTree, clustersTree, esd);
1828         if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1829       }
1830       if (fLoader[iDet]) {
1831         fLoader[iDet]->UnloadRecPoints();
1832       }
1833
1834       if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1835     }
1836   }
1837
1838   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1839     AliError(Form("the following detectors were not found: %s", 
1840                   detStr.Data()));
1841     if (fStopOnError) return kFALSE;
1842   }
1843   AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1844   eventNr++;
1845   return kTRUE;
1846 }
1847
1848 //_____________________________________________________________________________
1849 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1850 {
1851   // Reads the trigger decision which is
1852   // stored in Trigger.root file and fills
1853   // the corresponding esd entries
1854
1855   AliCodeTimerAuto("")
1856   
1857   AliInfo("Filling trigger information into the ESD");
1858
1859   AliCentralTrigger *aCTP = NULL;
1860
1861   if (fRawReader) {
1862     AliCTPRawStream input(fRawReader);
1863     if (!input.Next()) {
1864       AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1865       return kFALSE;
1866     }
1867     esd->SetTriggerMask(input.GetClassMask());
1868     esd->SetTriggerCluster(input.GetClusterMask());
1869
1870     aCTP = new AliCentralTrigger();
1871     TString configstr("");
1872     if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1873       AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1874       delete aCTP;
1875       return kFALSE;
1876     }
1877   }
1878   else {
1879     AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1880     if (runloader) {
1881       if (!runloader->LoadTrigger()) {
1882         aCTP = runloader->GetTrigger();
1883         esd->SetTriggerMask(aCTP->GetClassMask());
1884         esd->SetTriggerCluster(aCTP->GetClusterMask());
1885       }
1886       else {
1887         AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1888         return kFALSE;
1889       }
1890     }
1891     else {
1892       AliError("No run loader is available! The trigger information is not stored in the ESD !");
1893       return kFALSE;
1894     }
1895   }
1896
1897   // Now fill the trigger class names into AliESDRun object
1898   AliTriggerConfiguration *config = aCTP->GetConfiguration();
1899   if (!config) {
1900     AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1901     if (fRawReader) delete aCTP;
1902     return kFALSE;
1903   }
1904
1905   const TObjArray& classesArray = config->GetClasses();
1906   Int_t nclasses = classesArray.GetEntriesFast();
1907   for( Int_t j=0; j<nclasses; j++ ) {
1908     AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1909     Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1910     esd->SetTriggerClass(trclass->GetName(),trindex);
1911   }
1912
1913   if (fRawReader) delete aCTP;
1914   return kTRUE;
1915 }
1916
1917
1918
1919
1920
1921 //_____________________________________________________________________________
1922 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1923 {
1924   // 
1925   // Filling information from RawReader Header
1926   // 
1927
1928   AliInfo("Filling information from RawReader Header");
1929   esd->SetBunchCrossNumber(0);
1930   esd->SetOrbitNumber(0);
1931   esd->SetPeriodNumber(0);
1932   esd->SetTimeStamp(0);
1933   esd->SetEventType(0);
1934   const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1935   if (eventHeader){
1936
1937     const UInt_t *id = eventHeader->GetP("Id");
1938     esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1939     esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1940     esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1941
1942     esd->SetTimeStamp((eventHeader->Get("Timestamp")));  
1943     esd->SetEventType((eventHeader->Get("Type")));
1944   }
1945
1946   return kTRUE;
1947 }
1948
1949
1950 //_____________________________________________________________________________
1951 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1952 {
1953 // check whether detName is contained in detectors
1954 // if yes, it is removed from detectors
1955
1956   // check if all detectors are selected
1957   if ((detectors.CompareTo("ALL") == 0) ||
1958       detectors.BeginsWith("ALL ") ||
1959       detectors.EndsWith(" ALL") ||
1960       detectors.Contains(" ALL ")) {
1961     detectors = "ALL";
1962     return kTRUE;
1963   }
1964
1965   // search for the given detector
1966   Bool_t result = kFALSE;
1967   if ((detectors.CompareTo(detName) == 0) ||
1968       detectors.BeginsWith(detName+" ") ||
1969       detectors.EndsWith(" "+detName) ||
1970       detectors.Contains(" "+detName+" ")) {
1971     detectors.ReplaceAll(detName, "");
1972     result = kTRUE;
1973   }
1974
1975   // clean up the detectors string
1976   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1977   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1978   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1979
1980   return result;
1981 }
1982
1983 //_____________________________________________________________________________
1984 Bool_t AliReconstruction::InitRunLoader()
1985 {
1986 // get or create the run loader
1987
1988   if (gAlice) delete gAlice;
1989   gAlice = NULL;
1990
1991   if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1992     // load all base libraries to get the loader classes
1993     TString libs = gSystem->GetLibraries();
1994     for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1995       TString detName = fgkDetectorName[iDet];
1996       if (detName == "HLT") continue;
1997       if (libs.Contains("lib" + detName + "base.so")) continue;
1998       gSystem->Load("lib" + detName + "base.so");
1999     }
2000     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2001     if (!fRunLoader) {
2002       AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2003       CleanUp();
2004       return kFALSE;
2005     }
2006     fRunLoader->CdGAFile();
2007     if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
2008       if (fRunLoader->LoadgAlice() == 0) {
2009         gAlice = fRunLoader->GetAliRun();
2010         AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
2011       }
2012     }
2013     if (!gAlice && !fRawReader) {
2014       AliError(Form("no gAlice object found in file %s",
2015                     fGAliceFileName.Data()));
2016       CleanUp();
2017       return kFALSE;
2018     }
2019
2020     //PH This is a temporary fix to give access to the kinematics
2021     //PH that is needed for the labels of ITS clusters
2022     fRunLoader->LoadHeader();
2023     fRunLoader->LoadKinematics();
2024
2025   } else {               // galice.root does not exist
2026     if (!fRawReader) {
2027       AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2028       CleanUp();
2029       return kFALSE;
2030     }
2031     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2032                                     AliConfig::GetDefaultEventFolderName(),
2033                                     "recreate");
2034     if (!fRunLoader) {
2035       AliError(Form("could not create run loader in file %s", 
2036                     fGAliceFileName.Data()));
2037       CleanUp();
2038       return kFALSE;
2039     }
2040     fIsNewRunLoader = kTRUE;
2041     fRunLoader->MakeTree("E");
2042
2043     if (fNumberOfEventsPerFile > 0)
2044       fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2045     else
2046       fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2047   }
2048
2049   return kTRUE;
2050 }
2051
2052 //_____________________________________________________________________________
2053 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2054 {
2055 // get the reconstructor object and the loader for a detector
2056
2057   if (fReconstructor[iDet]) return fReconstructor[iDet];
2058
2059   // load the reconstructor object
2060   TPluginManager* pluginManager = gROOT->GetPluginManager();
2061   TString detName = fgkDetectorName[iDet];
2062   TString recName = "Ali" + detName + "Reconstructor";
2063
2064   if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2065
2066   AliReconstructor* reconstructor = NULL;
2067   // first check if a plugin is defined for the reconstructor
2068   TPluginHandler* pluginHandler = 
2069     pluginManager->FindHandler("AliReconstructor", detName);
2070   // if not, add a plugin for it
2071   if (!pluginHandler) {
2072     AliDebug(1, Form("defining plugin for %s", recName.Data()));
2073     TString libs = gSystem->GetLibraries();
2074     if (libs.Contains("lib" + detName + "base.so") ||
2075         (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2076       pluginManager->AddHandler("AliReconstructor", detName, 
2077                                 recName, detName + "rec", recName + "()");
2078     } else {
2079       pluginManager->AddHandler("AliReconstructor", detName, 
2080                                 recName, detName, recName + "()");
2081     }
2082     pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2083   }
2084   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2085     reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2086   }
2087   if (reconstructor) {
2088     TObject* obj = fOptions.FindObject(detName.Data());
2089     if (obj) reconstructor->SetOption(obj->GetTitle());
2090     reconstructor->Init();
2091     fReconstructor[iDet] = reconstructor;
2092   }
2093
2094   // get or create the loader
2095   if (detName != "HLT") {
2096     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2097     if (!fLoader[iDet]) {
2098       AliConfig::Instance()
2099         ->CreateDetectorFolders(fRunLoader->GetEventFolder(), 
2100                                 detName, detName);
2101       // first check if a plugin is defined for the loader
2102       pluginHandler = 
2103         pluginManager->FindHandler("AliLoader", detName);
2104       // if not, add a plugin for it
2105       if (!pluginHandler) {
2106         TString loaderName = "Ali" + detName + "Loader";
2107         AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2108         pluginManager->AddHandler("AliLoader", detName, 
2109                                   loaderName, detName + "base", 
2110                                   loaderName + "(const char*, TFolder*)");
2111         pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2112       }
2113       if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2114         fLoader[iDet] = 
2115           (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(), 
2116                                                  fRunLoader->GetEventFolder());
2117       }
2118       if (!fLoader[iDet]) {   // use default loader
2119         fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2120       }
2121       if (!fLoader[iDet]) {
2122         AliWarning(Form("couldn't get loader for %s", detName.Data()));
2123         if (fStopOnError) return NULL;
2124       } else {
2125         fRunLoader->AddLoader(fLoader[iDet]);
2126         fRunLoader->CdGAFile();
2127         if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2128         fRunLoader->Write(0, TObject::kOverwrite);
2129       }
2130     }
2131   }
2132       
2133   return reconstructor;
2134 }
2135
2136 //_____________________________________________________________________________
2137 Bool_t AliReconstruction::CreateVertexer()
2138 {
2139 // create the vertexer
2140
2141   fVertexer = NULL;
2142   AliReconstructor* itsReconstructor = GetReconstructor(0);
2143   if (itsReconstructor) {
2144     fVertexer = itsReconstructor->CreateVertexer();
2145   }
2146   if (!fVertexer) {
2147     AliWarning("couldn't create a vertexer for ITS");
2148     if (fStopOnError) return kFALSE;
2149   }
2150
2151   return kTRUE;
2152 }
2153
2154 //_____________________________________________________________________________
2155 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2156 {
2157 // create the trackers
2158
2159   TString detStr = detectors;
2160   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2161     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2162     AliReconstructor* reconstructor = GetReconstructor(iDet);
2163     if (!reconstructor) continue;
2164     TString detName = fgkDetectorName[iDet];
2165     if (detName == "HLT") {
2166       fRunHLTTracking = kTRUE;
2167       continue;
2168     }
2169     if (detName == "MUON") {
2170       fRunMuonTracking = kTRUE;
2171       continue;
2172     }
2173
2174
2175     fTracker[iDet] = reconstructor->CreateTracker();
2176     if (!fTracker[iDet] && (iDet < 7)) {
2177       AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2178       if (fStopOnError) return kFALSE;
2179     }
2180     AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2181   }
2182
2183   return kTRUE;
2184 }
2185
2186 //_____________________________________________________________________________
2187 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2188 {
2189 // delete trackers and the run loader and close and delete the file
2190
2191   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2192     delete fReconstructor[iDet];
2193     fReconstructor[iDet] = NULL;
2194     fLoader[iDet] = NULL;
2195     delete fTracker[iDet];
2196     fTracker[iDet] = NULL;
2197 //    delete fQADataMaker[iDet];
2198 //    fQADataMaker[iDet] = NULL;
2199   }
2200   delete fVertexer;
2201   fVertexer = NULL;
2202
2203   if (ftVertexer) delete ftVertexer;
2204   ftVertexer = NULL;
2205   
2206   if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2207         delete fDiamondProfile;
2208         fDiamondProfile = NULL;
2209         delete fDiamondProfileTPC;
2210         fDiamondProfileTPC = NULL;
2211         delete fGRPData;
2212         fGRPData = NULL;
2213   }
2214
2215
2216   delete fRunLoader;
2217   fRunLoader = NULL;
2218   delete fRawReader;
2219   fRawReader = NULL;
2220   if (fParentRawReader) delete fParentRawReader;
2221   fParentRawReader=NULL;
2222
2223   if (file) {
2224     file->Close();
2225     delete file;
2226   }
2227
2228   if (fileOld) {
2229     fileOld->Close();
2230     delete fileOld;
2231     gSystem->Unlink("AliESDs.old.root");
2232   }
2233
2234 }
2235
2236 //_____________________________________________________________________________
2237
2238 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2239 {
2240 // read the ESD event from a file
2241
2242   if (!esd) return kFALSE;
2243   char fileName[256];
2244   sprintf(fileName, "ESD_%d.%d_%s.root", 
2245           esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2246   if (gSystem->AccessPathName(fileName)) return kFALSE;
2247
2248   AliInfo(Form("reading ESD from file %s", fileName));
2249   AliDebug(1, Form("reading ESD from file %s", fileName));
2250   TFile* file = TFile::Open(fileName);
2251   if (!file || !file->IsOpen()) {
2252     AliError(Form("opening %s failed", fileName));
2253     delete file;
2254     return kFALSE;
2255   }
2256
2257   gROOT->cd();
2258   delete esd;
2259   esd = (AliESDEvent*) file->Get("ESD");
2260   file->Close();
2261   delete file;
2262   return kTRUE;
2263
2264 }
2265
2266
2267
2268 //_____________________________________________________________________________
2269 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2270 {
2271 // write the ESD event to a file
2272
2273   if (!esd) return;
2274   char fileName[256];
2275   sprintf(fileName, "ESD_%d.%d_%s.root", 
2276           esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2277
2278   AliDebug(1, Form("writing ESD to file %s", fileName));
2279   TFile* file = TFile::Open(fileName, "recreate");
2280   if (!file || !file->IsOpen()) {
2281     AliError(Form("opening %s failed", fileName));
2282   } else {
2283     esd->Write("ESD");
2284     file->Close();
2285   }
2286   delete file;
2287 }
2288
2289
2290 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2291 {
2292   // Write space-points which are then used in the alignment procedures
2293   // For the moment only ITS, TRD and TPC
2294
2295   // Load TOF clusters
2296   if (fTracker[3]){
2297     fLoader[3]->LoadRecPoints("read");
2298     TTree* tree = fLoader[3]->TreeR();
2299     if (!tree) {
2300       AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2301       return;
2302     }
2303     fTracker[3]->LoadClusters(tree);
2304   }
2305   Int_t ntracks = esd->GetNumberOfTracks();
2306   for (Int_t itrack = 0; itrack < ntracks; itrack++)
2307     {
2308       AliESDtrack *track = esd->GetTrack(itrack);
2309       Int_t nsp = 0;
2310       Int_t idx[200];
2311       for (Int_t iDet = 3; iDet >= 0; iDet--)
2312         nsp += track->GetNcls(iDet);
2313       if (nsp) {
2314         AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2315         track->SetTrackPointArray(sp);
2316         Int_t isptrack = 0;
2317         for (Int_t iDet = 3; iDet >= 0; iDet--) {
2318           AliTracker *tracker = fTracker[iDet];
2319           if (!tracker) continue;
2320           Int_t nspdet = track->GetNcls(iDet);
2321           if (nspdet <= 0) continue;
2322           track->GetClusters(iDet,idx);
2323           AliTrackPoint p;
2324           Int_t isp = 0;
2325           Int_t isp2 = 0;
2326           while (isp2 < nspdet) {
2327             Bool_t isvalid;
2328             TString dets = fgkDetectorName[iDet];
2329             if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2330             fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2331             fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2332             fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2333               isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2334             } else {
2335               isvalid = tracker->GetTrackPoint(idx[isp2],p); 
2336             } 
2337             isp2++;
2338             const Int_t kNTPCmax = 159;
2339             if (iDet==1 && isp2>kNTPCmax) break;   // to be fixed
2340             if (!isvalid) continue;
2341             sp->AddPoint(isptrack,&p); isptrack++; isp++;
2342           }
2343         }       
2344       }
2345     }
2346   if (fTracker[3]){
2347     fTracker[3]->UnloadClusters();
2348     fLoader[3]->UnloadRecPoints();
2349   }
2350 }
2351
2352 //_____________________________________________________________________________
2353 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2354 {
2355   // The method reads the raw-data error log
2356   // accumulated within the rawReader.
2357   // It extracts the raw-data errors related to
2358   // the current event and stores them into
2359   // a TClonesArray inside the esd object.
2360
2361   if (!fRawReader) return;
2362
2363   for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2364
2365     AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2366     if (!log) continue;
2367     if (iEvent != log->GetEventNumber()) continue;
2368
2369     esd->AddRawDataErrorLog(log);
2370   }
2371
2372 }
2373
2374 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2375   // Dump a file content into a char in TNamed
2376   ifstream in;
2377   in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2378   Int_t kBytes = (Int_t)in.tellg();
2379   printf("Size: %d \n",kBytes);
2380   TNamed *fn = 0;
2381   if(in.good()){
2382     char* memblock = new char [kBytes];
2383     in.seekg (0, ios::beg);
2384     in.read (memblock, kBytes);
2385     in.close();
2386     TString fData(memblock,kBytes);
2387     fn = new TNamed(pName,fData);
2388     printf("fData Size: %d \n",fData.Sizeof());
2389     printf("pName Size: %d \n",pName.Sizeof());
2390     printf("fn    Size: %d \n",fn->Sizeof());
2391     delete[] memblock;
2392   }
2393   else{
2394     AliInfo(Form("Could not Open %s\n",fPath.Data()));
2395   }
2396
2397   return fn;
2398 }
2399
2400 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2401   // This is not really needed in AliReconstruction at the moment
2402   // but can serve as a template
2403
2404   TList *fList = fTree->GetUserInfo();
2405   TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2406   printf("fn Size: %d \n",fn->Sizeof());
2407
2408   TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2409   const char* cdata = fn->GetTitle();
2410   printf("fTmp Size %d\n",fTmp.Sizeof());
2411
2412   int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2413   printf("calculated size %d\n",size);
2414   ofstream out(pName.Data(),ios::out | ios::binary);
2415   out.write(cdata,size);
2416   out.close();
2417
2418 }
2419   
2420 //_____________________________________________________________________________
2421 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2422 {
2423  // get the quality assurance data maker object and the loader for a detector
2424
2425   if (fQADataMaker[iDet]) 
2426     return fQADataMaker[iDet];
2427
2428   AliQADataMakerRec * qadm = NULL;
2429   if (iDet == fgkNDetectors) { //Global QA
2430      qadm = new AliGlobalQADataMaker();
2431      fQADataMaker[iDet] = qadm;
2432      return qadm;
2433   }
2434
2435   // load the QA data maker object
2436   TPluginManager* pluginManager = gROOT->GetPluginManager();
2437   TString detName = fgkDetectorName[iDet];
2438   TString qadmName = "Ali" + detName + "QADataMakerRec";
2439   if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) 
2440     return NULL;
2441
2442   // first check if a plugin is defined for the quality assurance data maker
2443   TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2444   // if not, add a plugin for it
2445   if (!pluginHandler) {
2446     AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2447     TString libs = gSystem->GetLibraries();
2448     if (libs.Contains("lib" + detName + "base.so") ||
2449         (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2450       pluginManager->AddHandler("AliQADataMakerRec", detName, 
2451                                 qadmName, detName + "qadm", qadmName + "()");
2452     } else {
2453       pluginManager->AddHandler("AliQADataMakerRec", detName, 
2454                                 qadmName, detName, qadmName + "()");
2455     }
2456     pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2457   }
2458   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2459     qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2460   }
2461
2462   fQADataMaker[iDet] = qadm;
2463
2464   return qadm;
2465 }
2466
2467 //_____________________________________________________________________________
2468 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
2469 {
2470   // run the Quality Assurance data producer
2471
2472   AliCodeTimerAuto("")
2473   TString detStr = detectors;
2474   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2475    if (!IsSelected(fgkDetectorName[iDet], detStr)) 
2476      continue;
2477    AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2478    if (!qadm) 
2479      continue;
2480    AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2481    AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2482     
2483    qadm->Exec(AliQA::kESDS, esd) ; 
2484    qadm->Increment() ; 
2485
2486    AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2487  }
2488  if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2489    AliError(Form("the following detectors were not found: %s",
2490                  detStr.Data()));
2491    if (fStopOnError) 
2492      return kFALSE;
2493  }
2494  
2495  return kTRUE;
2496   
2497 }
2498
2499 //_____________________________________________________________________________
2500 void AliReconstruction::CheckQA()
2501 {
2502 // check the QA of SIM for this run and remove the detectors 
2503 // with status Fatal
2504   
2505         TString newRunLocalReconstruction ; 
2506         TString newRunTracking ;
2507         TString newFillESD ;
2508          
2509         for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2510                 TString detName(AliQA::GetDetName(iDet)) ;
2511                 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ; 
2512                 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2513                                 AliInfo(Form("QA status for %s in Hits and/or SDIGITS  and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2514                 } else {
2515                         if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) || 
2516                                         fRunLocalReconstruction.Contains("ALL") )  {
2517                                 newRunLocalReconstruction += detName ; 
2518                                 newRunLocalReconstruction += " " ;                      
2519                         }
2520                         if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) || 
2521                                         fRunTracking.Contains("ALL") )  {
2522                                 newRunTracking += detName ; 
2523                                 newRunTracking += " " ;                         
2524                         }
2525                         if ( fFillESD.Contains(AliQA::GetDetName(iDet)) || 
2526                                         fFillESD.Contains("ALL") )  {
2527                                 newFillESD += detName ; 
2528                                 newFillESD += " " ;                     
2529                         }
2530                 }
2531         }
2532         fRunLocalReconstruction = newRunLocalReconstruction ; 
2533         fRunTracking            = newRunTracking ; 
2534         fFillESD                = newFillESD ; 
2535 }
2536
2537 //_____________________________________________________________________________
2538 Int_t AliReconstruction::GetDetIndex(const char* detector)
2539 {
2540   // return the detector index corresponding to detector
2541   Int_t index = -1 ; 
2542   for (index = 0; index < fgkNDetectors ; index++) {
2543     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2544         break ; 
2545   }     
2546   return index ; 
2547 }
2548 //_____________________________________________________________________________
2549 Bool_t AliReconstruction::FinishPlaneEff() {
2550  //
2551  // Here execute all the necessary operationis, at the end of the tracking phase,
2552  // in case that evaluation of PlaneEfficiencies was required for some detector. 
2553  // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated. 
2554  //
2555  // This Preliminary version works only FOR ITS !!!!!
2556  // other detectors (TOF,TRD, etc. have to develop their specific codes)
2557  //
2558  //  Input: none
2559  //  Return: kTRUE if all operations have been done properly, kFALSE otherwise 
2560  //
2561  Bool_t ret=kFALSE;
2562  //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2563  for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS  
2564    //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2565    if(fTracker[iDet]) {
2566       AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff(); 
2567       ret=planeeff->WriteIntoCDB();
2568       if(planeeff->GetCreateHistos()) {
2569         TString name="PlaneEffHisto";
2570         name+=fgkDetectorName[iDet];
2571         name+=".root";
2572         ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2573       }
2574    }
2575  }
2576  return ret;
2577 }
2578 //_____________________________________________________________________________
2579 Bool_t AliReconstruction::InitPlaneEff() {
2580 //
2581  // Here execute all the necessary operations, before of the tracking phase,
2582  // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2583  // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency 
2584  // which should be updated/recalculated.
2585  //
2586  // This Preliminary version will work only FOR ITS !!!!!
2587  // other detectors (TOF,TRD, etc. have to develop their specific codes)
2588  //
2589  //  Input: none
2590  //  Return: kTRUE if all operations have been done properly, kFALSE otherwise
2591  //
2592  AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2593  return kTRUE;
2594 }
2595
2596 //_____________________________________________________________________________
2597 Bool_t AliReconstruction::InitAliEVE()
2598 {
2599   // This method should be called only in case 
2600   // AliReconstruction is run
2601   // within the alieve environment.
2602   // It will initialize AliEVE in a way
2603   // so that it can visualize event processed
2604   // by AliReconstruction.
2605   // The return flag shows whenever the
2606   // AliEVE initialization was successful or not.
2607
2608   TString macroStr;
2609   macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2610   AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2611   if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2612
2613   gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online()\");gEve->AddEvent(gAliEveEvent);};");
2614
2615   return kTRUE;
2616 }
2617   
2618 //_____________________________________________________________________________
2619 void AliReconstruction::RunAliEVE()
2620 {
2621   // Runs AliEVE visualisation of
2622   // the current event.
2623   // Should be executed only after
2624   // successful initialization of AliEVE.
2625
2626   AliInfo("Running AliEVE...");
2627   gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2628   gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2629   gSystem->Run();
2630 }