]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.cxx
Splitting of the QA maker into simulation and reconstruction dependent parts (Yves)
[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 // For debug purposes the method SetCheckPointLevel can be used. If the      //
105 // argument is greater than 0, files with ESD events will be written after   //
106 // selected steps of the reconstruction for each event:                      //
107 //   level 1: after tracking and after filling of ESD (final)                //
108 //   level 2: in addition after each tracking step                           //
109 //   level 3: in addition after the filling of ESD for each detector         //
110 // If a final check point file exists for an event, this event will be       //
111 // skipped in the reconstruction. The tracking and the filling of ESD for    //
112 // a detector will be skipped as well, if the corresponding check point      //
113 // file exists. The ESD event will then be loaded from the file instead.     //
114 //                                                                           //
115 ///////////////////////////////////////////////////////////////////////////////
116
117 #include <TArrayF.h>
118 #include <TFile.h>
119 #include <TList.h>
120 #include <TSystem.h>
121 #include <TROOT.h>
122 #include <TPluginManager.h>
123 #include <TGeoManager.h>
124 #include <TLorentzVector.h>
125 #include <TArrayS.h>
126 #include <TArrayD.h>
127
128 #include "AliReconstruction.h"
129 #include "AliCodeTimer.h"
130 #include "AliReconstructor.h"
131 #include "AliLog.h"
132 #include "AliRunLoader.h"
133 #include "AliRun.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderDate.h"
136 #include "AliRawReaderRoot.h"
137 #include "AliRawEventHeaderBase.h"
138 #include "AliESDEvent.h"
139 #include "AliESDMuonTrack.h"
140 #include "AliESDfriend.h"
141 #include "AliESDVertex.h"
142 #include "AliESDcascade.h"
143 #include "AliESDkink.h"
144 #include "AliESDtrack.h"
145 #include "AliESDCaloCluster.h"
146 #include "AliESDCaloCells.h"
147 #include "AliMultiplicity.h"
148 #include "AliTracker.h"
149 #include "AliVertexer.h"
150 #include "AliVertexerTracks.h"
151 #include "AliV0vertexer.h"
152 #include "AliCascadeVertexer.h"
153 #include "AliHeader.h"
154 #include "AliGenEventHeader.h"
155 #include "AliPID.h"
156 #include "AliESDpid.h"
157 #include "AliESDtrack.h"
158 #include "AliESDPmdTrack.h"
159
160 #include "AliESDTagCreator.h"
161 #include "AliAODTagCreator.h"
162
163 #include "AliGeomManager.h"
164 #include "AliTrackPointArray.h"
165 #include "AliCDBManager.h"
166 #include "AliCDBStorage.h"
167 #include "AliCDBEntry.h"
168 #include "AliAlignObj.h"
169
170 #include "AliCentralTrigger.h"
171 #include "AliCTPRawStream.h"
172
173 #include "AliAODEvent.h"
174 #include "AliAODHeader.h"
175 #include "AliAODTrack.h"
176 #include "AliAODVertex.h"
177 #include "AliAODv0.h"
178 #include "AliAODJet.h"
179 #include "AliAODCaloCells.h"
180 #include "AliAODCaloCluster.h"
181 #include "AliAODPmdCluster.h"
182 #include "AliAODFmdCluster.h"
183 #include "AliAODTracklets.h"
184
185 #include "AliQADataMakerRec.h" 
186 #include "AliGlobalQADataMaker.h" 
187 #include "AliQA.h"
188 #include "AliQADataMakerSteer.h"
189
190 #include "AliSysInfo.h" // memory snapshots
191
192
193 ClassImp(AliReconstruction)
194
195
196 //_____________________________________________________________________________
197 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
198
199 //_____________________________________________________________________________
200 AliReconstruction::AliReconstruction(const char* gAliceFilename,
201                                      const char* name, const char* title) :
202   TNamed(name, title),
203
204   fUniformField(kTRUE),
205   fRunVertexFinder(kTRUE),
206   fRunHLTTracking(kFALSE),
207   fRunMuonTracking(kFALSE),
208   fRunV0Finder(kTRUE),
209   fRunCascadeFinder(kTRUE),
210   fStopOnError(kFALSE),
211   fWriteAlignmentData(kFALSE),
212   fWriteESDfriend(kFALSE),
213   fWriteAOD(kFALSE),
214   fFillTriggerESD(kTRUE),
215
216   fCleanESD(kTRUE),
217   fV0DCAmax(3.),
218   fV0CsPmin(0.),
219   fDmax(50.),
220   fZmax(50.),
221
222   fRunLocalReconstruction("ALL"),
223   fRunTracking("ALL"),
224   fFillESD("ALL"),
225   fUseTrackingErrorsForAlignment(""),
226   fGAliceFileName(gAliceFilename),
227   fInput(""),
228   fEquipIdMap(""),
229   fFirstEvent(0),
230   fLastEvent(-1),
231   fNumberOfEventsPerFile(1),
232   fCheckPointLevel(0),
233   fOptions(),
234   fLoadAlignFromCDB(kTRUE),
235   fLoadAlignData("ALL"),
236   fESDPar(""),
237
238   fRunLoader(NULL),
239   fRawReader(NULL),
240
241   fVertexer(NULL),
242   fDiamondProfile(NULL),
243   fMeanVertexConstraint(kTRUE),
244
245   fGRPList(NULL),
246
247   fAlignObjArray(NULL),
248   fCDBUri(),
249   fSpecCDBUri(), 
250   fInitCDBCalled(kFALSE),
251   fSetRunNumberFromDataCalled(kFALSE),
252   fRunQA(kTRUE),  
253   fRunGlobalQA(kFALSE),
254   fInLoopQA(kFALSE)
255 {
256 // create reconstruction object with default parameters
257   
258   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
259     fReconstructor[iDet] = NULL;
260     fLoader[iDet] = NULL;
261     fTracker[iDet] = NULL;
262     fQADataMaker[iDet] = NULL;
263         fQACycles[iDet] = 999999;       
264   }
265   fQADataMaker[fgkNDetectors]=NULL;  //Global QA
266   AliPID pid;
267 }
268
269 //_____________________________________________________________________________
270 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
271   TNamed(rec),
272
273   fUniformField(rec.fUniformField),
274   fRunVertexFinder(rec.fRunVertexFinder),
275   fRunHLTTracking(rec.fRunHLTTracking),
276   fRunMuonTracking(rec.fRunMuonTracking),
277   fRunV0Finder(rec.fRunV0Finder),
278   fRunCascadeFinder(rec.fRunCascadeFinder),
279   fStopOnError(rec.fStopOnError),
280   fWriteAlignmentData(rec.fWriteAlignmentData),
281   fWriteESDfriend(rec.fWriteESDfriend),
282   fWriteAOD(rec.fWriteAOD),
283   fFillTriggerESD(rec.fFillTriggerESD),
284
285   fCleanESD(rec.fCleanESD),
286   fV0DCAmax(rec.fV0DCAmax),
287   fV0CsPmin(rec.fV0CsPmin),
288   fDmax(rec.fDmax),
289   fZmax(rec.fZmax),
290
291   fRunLocalReconstruction(rec.fRunLocalReconstruction),
292   fRunTracking(rec.fRunTracking),
293   fFillESD(rec.fFillESD),
294   fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
295   fGAliceFileName(rec.fGAliceFileName),
296   fInput(rec.fInput),
297   fEquipIdMap(rec.fEquipIdMap),
298   fFirstEvent(rec.fFirstEvent),
299   fLastEvent(rec.fLastEvent),
300   fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
301   fCheckPointLevel(0),
302   fOptions(),
303   fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
304   fLoadAlignData(rec.fLoadAlignData),
305   fESDPar(rec.fESDPar),
306
307   fRunLoader(NULL),
308   fRawReader(NULL),
309
310   fVertexer(NULL),
311   fDiamondProfile(NULL),
312   fMeanVertexConstraint(rec.fMeanVertexConstraint),
313
314   fGRPList(NULL),
315
316   fAlignObjArray(rec.fAlignObjArray),
317   fCDBUri(rec.fCDBUri),
318   fSpecCDBUri(), 
319   fInitCDBCalled(rec.fInitCDBCalled),
320   fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
321   fRunQA(rec.fRunQA),  
322   fRunGlobalQA(rec.fRunGlobalQA),
323   fInLoopQA(rec.fInLoopQA)
324 {
325 // copy constructor
326
327   for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
328     if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
329   }
330   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
331     fReconstructor[iDet] = NULL;
332     fLoader[iDet] = NULL;
333     fTracker[iDet] = NULL;
334     fQADataMaker[iDet] = NULL;
335         fQACycles[iDet] = rec.fQACycles[iDet];  
336   }
337   fQADataMaker[fgkNDetectors]=NULL;  //Global QA
338   for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
339     if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
340   }
341 }
342
343 //_____________________________________________________________________________
344 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
345 {
346 // assignment operator
347
348   this->~AliReconstruction();
349   new(this) AliReconstruction(rec);
350   return *this;
351 }
352
353 //_____________________________________________________________________________
354 AliReconstruction::~AliReconstruction()
355 {
356 // clean up
357
358   CleanUp();
359   fOptions.Delete();
360   fSpecCDBUri.Delete();
361
362   AliCodeTimer::Instance()->Print();
363 }
364
365 //_____________________________________________________________________________
366 void AliReconstruction::InitCDB()
367 {
368 // activate a default CDB storage
369 // First check if we have any CDB storage set, because it is used 
370 // to retrieve the calibration and alignment constants
371
372   if (fInitCDBCalled) return;
373   fInitCDBCalled = kTRUE;
374
375   AliCDBManager* man = AliCDBManager::Instance();
376   if (man->IsDefaultStorageSet())
377   {
378     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
379     AliWarning("Default CDB storage has been already set !");
380     AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
381     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
382     fCDBUri = man->GetDefaultStorage()->GetURI();
383   }
384   else {
385     if (fCDBUri.Length() > 0) 
386     {
387         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
388         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
389         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
390     } else {
391         fCDBUri="local://$ALICE_ROOT";
392         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
393         AliWarning("Default CDB storage not yet set !!!!");
394         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
395         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
396                 
397     }
398     man->SetDefaultStorage(fCDBUri);
399   }
400
401   // Now activate the detector specific CDB storage locations
402   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
403     TObject* obj = fSpecCDBUri[i];
404     if (!obj) continue;
405     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
406     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
407     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
408     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
409   }
410   
411 }
412
413 //_____________________________________________________________________________
414 void AliReconstruction::SetDefaultStorage(const char* uri) {
415 // Store the desired default CDB storage location
416 // Activate it later within the Run() method
417
418   fCDBUri = uri;
419
420 }
421
422 //_____________________________________________________________________________
423 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
424 // Store a detector-specific CDB storage location
425 // Activate it later within the Run() method
426
427   AliCDBPath aPath(calibType);
428   if(!aPath.IsValid()){
429         // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
430         for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
431                 if(!strcmp(calibType, fgkDetectorName[iDet])) {
432                         aPath.SetPath(Form("%s/*", calibType));
433                         AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
434                         break;
435                 }
436         }
437         if(!aPath.IsValid()){
438                 AliError(Form("Not a valid path or detector: %s", calibType));
439                 return;
440         }
441   }
442
443 //  // check that calibType refers to a "valid" detector name
444 //  Bool_t isDetector = kFALSE;
445 //  for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
446 //    TString detName = fgkDetectorName[iDet];
447 //    if(aPath.GetLevel0() == detName) {
448 //      isDetector = kTRUE;
449 //      break;
450 //    }
451 //  }
452 //
453 //  if(!isDetector) {
454 //      AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
455 //      return;
456 //  }
457
458   TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
459   if (obj) fSpecCDBUri.Remove(obj);
460   fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
461
462 }
463
464 //_____________________________________________________________________________
465 Bool_t AliReconstruction::SetRunNumberFromData()
466 {
467   // The method is called in Run() in order
468   // to set a correct run number.
469   // In case of raw data reconstruction the
470   // run number is taken from the raw data header
471
472   if (fSetRunNumberFromDataCalled) return kTRUE;
473   fSetRunNumberFromDataCalled = kTRUE;
474   
475   AliCDBManager* man = AliCDBManager::Instance();
476   
477   if(man->GetRun() > 0) {
478         AliWarning("Run number is taken from event header! Ignoring settings in AliCDBManager!");
479   } 
480   
481   if (!fRunLoader) {
482       AliError("No run loader is found !"); 
483       return kFALSE;
484     }
485     // read run number from gAlice
486     if(fRunLoader->GetAliRun())
487       AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
488     else {
489       if(fRawReader) {
490         if(fRawReader->NextEvent()) {
491           AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
492           fRawReader->RewindEvents();
493         }
494         else {
495           AliError("No raw-data events found !");
496           return kFALSE;
497         }
498       }
499       else {
500         AliError("Neither gAlice nor RawReader objects are found !");
501         return kFALSE;
502       }
503   }
504
505   man->Print();  
506   
507   return kTRUE;
508 }
509
510 //_____________________________________________________________________________
511 void AliReconstruction::SetCDBLock() {
512   // Set CDB lock: from now on it is forbidden to reset the run number
513   // or the default storage or to activate any further storage!
514   
515   AliCDBManager::Instance()->SetLock(1);
516 }
517
518 //_____________________________________________________________________________
519 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
520 {
521   // Read the alignment objects from CDB.
522   // Each detector is supposed to have the
523   // alignment objects in DET/Align/Data CDB path.
524   // All the detector objects are then collected,
525   // sorted by geometry level (starting from ALIC) and
526   // then applied to the TGeo geometry.
527   // Finally an overlaps check is performed.
528
529   // Load alignment data from CDB and fill fAlignObjArray 
530   if(fLoadAlignFromCDB){
531         
532     TString detStr = detectors;
533     TString loadAlObjsListOfDets = "";
534     
535     for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
536       if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
537       loadAlObjsListOfDets += fgkDetectorName[iDet];
538       loadAlObjsListOfDets += " ";
539     } // end loop over detectors
540     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
541     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
542   }else{
543     // Check if the array with alignment objects was
544     // provided by the user. If yes, apply the objects
545     // to the present TGeo geometry
546     if (fAlignObjArray) {
547       if (gGeoManager && gGeoManager->IsClosed()) {
548         if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
549           AliError("The misalignment of one or more volumes failed!"
550                    "Compare the list of simulated detectors and the list of detector alignment data!");
551           return kFALSE;
552         }
553       }
554       else {
555         AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
556         return kFALSE;
557       }
558     }
559   }
560   
561   delete fAlignObjArray; fAlignObjArray=0;
562
563   return kTRUE;
564 }
565
566 //_____________________________________________________________________________
567 void AliReconstruction::SetGAliceFile(const char* fileName)
568 {
569 // set the name of the galice file
570
571   fGAliceFileName = fileName;
572 }
573
574 //_____________________________________________________________________________
575 void AliReconstruction::SetOption(const char* detector, const char* option)
576 {
577 // set options for the reconstruction of a detector
578
579   TObject* obj = fOptions.FindObject(detector);
580   if (obj) fOptions.Remove(obj);
581   fOptions.Add(new TNamed(detector, option));
582 }
583
584
585 //_____________________________________________________________________________
586 Bool_t AliReconstruction::Run(const char* input)
587 {
588 // run the reconstruction
589
590   AliCodeTimerAuto("")
591   
592   // set the input
593   if (!input) input = fInput.Data();
594   TString fileName(input);
595   if (fileName.EndsWith("/")) {
596     fRawReader = new AliRawReaderFile(fileName);
597   } else if (fileName.EndsWith(".root")) {
598     fRawReader = new AliRawReaderRoot(fileName);
599   } else if (!fileName.IsNull()) {
600     fRawReader = new AliRawReaderDate(fileName);
601     fRawReader->SelectEvents(7);
602   }
603   if (!fEquipIdMap.IsNull() && fRawReader)
604     fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
605
606    AliSysInfo::AddStamp("Start");
607   // get the run loader
608   if (!InitRunLoader()) return kFALSE;
609    AliSysInfo::AddStamp("LoadLoader");
610
611   // Initialize the CDB storage
612   InitCDB();
613   
614   AliSysInfo::AddStamp("LoadCDB");
615
616   // Set run number in CDBManager (if it is not already set by the user)
617   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
618   
619   // Set CDB lock: from now on it is forbidden to reset the run number
620   // or the default storage or to activate any further storage!
621   SetCDBLock();
622   
623   // Import ideal TGeo geometry and apply misalignment
624   if (!gGeoManager) {
625     TString geom(gSystem->DirName(fGAliceFileName));
626     geom += "/geometry.root";
627     AliGeomManager::LoadGeometry(geom.Data());
628     if (!gGeoManager) if (fStopOnError) return kFALSE;
629   }
630
631   if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
632    AliSysInfo::AddStamp("LoadGeom");
633
634   //QA 
635   AliQADataMakerSteer qas ; 
636   if (fRunQA && fRawReader) qas.Run(fRunLocalReconstruction, fRawReader) ; 
637   // checking the QA of previous steps
638   //CheckQA() ; 
639  
640   /*
641   // local reconstruction
642   if (!fRunLocalReconstruction.IsNull()) {
643     if (!RunLocalReconstruction(fRunLocalReconstruction)) {
644       if (fStopOnError) {CleanUp(); return kFALSE;}
645     }
646   }
647   */
648
649   // get vertexer
650   if (fRunVertexFinder && !CreateVertexer()) {
651     if (fStopOnError) {
652       CleanUp(); 
653       return kFALSE;
654     }
655   }
656    AliSysInfo::AddStamp("Vertexer");
657
658   // get trackers
659   if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
660     if (fStopOnError) {
661       CleanUp(); 
662       return kFALSE;
663     }      
664   }
665    AliSysInfo::AddStamp("LoadTrackers");
666
667   // get the possibly already existing ESD file and tree
668   AliESDEvent* esd = new AliESDEvent(); AliESDEvent* hltesd = new AliESDEvent();
669   TFile* fileOld = NULL;
670   TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
671   if (!gSystem->AccessPathName("AliESDs.root")){
672     gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
673     fileOld = TFile::Open("AliESDs.old.root");
674     if (fileOld && fileOld->IsOpen()) {
675       treeOld = (TTree*) fileOld->Get("esdTree");
676       if (treeOld)esd->ReadFromTree(treeOld);
677       hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
678       if (hlttreeOld)   hltesd->ReadFromTree(hlttreeOld);
679     }
680   }
681
682   // create the ESD output file and tree
683   TFile* file = TFile::Open("AliESDs.root", "RECREATE");
684   file->SetCompressionLevel(2);
685   if (!file->IsOpen()) {
686     AliError("opening AliESDs.root failed");
687     if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}    
688   }
689
690   TTree* tree = new TTree("esdTree", "Tree with ESD objects");
691   esd = new AliESDEvent();
692   esd->CreateStdContent();
693   esd->WriteToTree(tree);
694
695   TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
696   hltesd = new AliESDEvent();
697   hltesd->CreateStdContent();
698   hltesd->WriteToTree(hlttree);
699
700   /* CKB Why?
701   delete esd; delete hltesd;
702   esd = NULL; hltesd = NULL;
703   */
704   // create the branch with ESD additions
705
706
707
708   AliESDfriend *esdf = 0; 
709   if (fWriteESDfriend) {
710     esdf = new AliESDfriend();
711     TBranch *br=tree->Branch("ESDfriend.","AliESDfriend", &esdf);
712     br->SetFile("AliESDfriends.root");
713     esd->AddObject(esdf);
714   }
715
716   
717   // Get the GRP CDB entry
718   AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
719         
720   if(entryGRP) {
721         fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());  
722   } else {
723         AliError("No GRP entry found in OCDB!");
724   }
725
726   // Get the diamond profile from OCDB
727   AliCDBEntry* entry = AliCDBManager::Instance()
728         ->Get("GRP/Calib/MeanVertex");
729         
730   if(entry) {
731         fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());  
732   } else {
733         AliError("No diamond profile found in OCDB!");
734   }
735
736   AliVertexerTracks tVertexer(AliTracker::GetBz());
737   if(fDiamondProfile && fMeanVertexConstraint) tVertexer.SetVtxStart(fDiamondProfile);
738
739   if (fRawReader) fRawReader->RewindEvents();
740
741   ProcInfo_t ProcInfo;
742   gSystem->GetProcInfo(&ProcInfo);
743   AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
744   
745
746   //Initialize the QA and start of cycle for out-of-cycle QA
747   if (fRunQA) {
748      TString detStr(fFillESD); 
749      for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
750         if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
751         AliQADataMakerRec *qadm = GetQADataMaker(iDet);  
752         if (!qadm) continue;
753         AliInfo(Form("Initializing the QA data maker for %s", 
754                fgkDetectorName[iDet]));
755         qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
756         qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
757         if (!fInLoopQA) {
758            qadm->StartOfCycle(AliQA::kRECPOINTS);
759            qadm->StartOfCycle(AliQA::kESDS,"same");
760         }
761      }
762      if (fRunGlobalQA) {
763         AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
764         AliInfo(Form("Initializing the global QA data maker"));
765         TObjArray *arr=
766           qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
767         AliTracker::SetResidualsArray(arr);
768         if (!fInLoopQA) {
769            qadm->StartOfCycle(AliQA::kRECPOINTS);
770         }
771      }
772   }
773
774
775   //******* The loop over events
776   for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
777     if (fRawReader) fRawReader->NextEvent();
778     if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
779       // copy old ESD to the new one
780       if (treeOld) {
781         esd->ReadFromTree(treeOld);
782         treeOld->GetEntry(iEvent);
783       }
784       tree->Fill();
785       if (hlttreeOld) {
786         esd->ReadFromTree(hlttreeOld);
787         hlttreeOld->GetEntry(iEvent);
788       }
789       hlttree->Fill();
790       continue;
791     }
792     
793     AliInfo(Form("processing event %d", iEvent));
794
795     //Start of cycle for the in-loop QA
796     if (fRunQA && fInLoopQA) {
797        TString detStr(fFillESD); 
798        for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
799           if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
800           AliQADataMakerRec *qadm = GetQADataMaker(iDet);  
801           if (!qadm) continue;
802           qadm->StartOfCycle(AliQA::kRECPOINTS);
803           qadm->StartOfCycle(AliQA::kESDS, "same") ;    
804        }
805        if (fRunGlobalQA) {
806           AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
807           qadm->StartOfCycle(AliQA::kRECPOINTS);
808        }
809     }
810
811     fRunLoader->GetEvent(iEvent);
812
813     char aFileName[256];
814     sprintf(aFileName, "ESD_%d.%d_final.root", 
815             fRunLoader->GetHeader()->GetRun(), 
816             fRunLoader->GetHeader()->GetEventNrInRun());
817     if (!gSystem->AccessPathName(aFileName)) continue;
818
819     // local signle event reconstruction
820     if (!fRunLocalReconstruction.IsNull()) {
821       if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
822         if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
823       }
824     }
825
826     esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
827     hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
828     esd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
829     hltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
830     
831     // Set magnetic field from the tracker
832     esd->SetMagneticField(AliTracker::GetBz());
833     hltesd->SetMagneticField(AliTracker::GetBz());
834
835     
836     
837     // Fill raw-data error log into the ESD
838     if (fRawReader) FillRawDataErrorLog(iEvent,esd);
839
840     // vertex finder
841     if (fRunVertexFinder) {
842       if (!ReadESD(esd, "vertex")) {
843         if (!RunVertexFinder(esd)) {
844           if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
845         }
846         if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
847       }
848     }
849
850     // HLT tracking
851     if (!fRunTracking.IsNull()) {
852       if (fRunHLTTracking) {
853         hltesd->SetVertex(esd->GetVertex());
854         if (!RunHLTTracking(hltesd)) {
855           if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
856         }
857       }
858     }
859
860     // Muon tracking
861     if (!fRunTracking.IsNull()) {
862       if (fRunMuonTracking) {
863         if (!RunMuonTracking(esd)) {
864           if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
865         }
866       }
867     }
868
869     // barrel tracking
870     if (!fRunTracking.IsNull()) {
871       if (!ReadESD(esd, "tracking")) {
872         if (!RunTracking(esd)) {
873           if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
874         }
875         if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
876       }
877     }
878
879     // fill ESD
880     if (!fFillESD.IsNull()) {
881       if (!FillESD(esd, fFillESD)) {
882         if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
883       }
884     }
885   
886     // fill Event header information from the RawEventHeader
887     if (fRawReader){FillRawEventHeaderESD(esd);}
888
889     // combined PID
890     AliESDpid::MakePID(esd);
891     if (fCheckPointLevel > 1) WriteESD(esd, "PID");
892
893     if (fFillTriggerESD) {
894       if (!ReadESD(esd, "trigger")) {
895         if (!FillTriggerESD(esd)) {
896           if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
897         }
898         if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
899       }
900     }
901
902     file->cd();
903
904     //Try to improve the reconstructed primary vertex position using the tracks
905     AliESDVertex *pvtx=0;
906     Bool_t dovertex=kTRUE;
907     TObject* obj = fOptions.FindObject("ITS");
908     if (obj) {
909       TString optITS = obj->GetTitle();
910       if (optITS.Contains("cosmics") || optITS.Contains("COSMICS")) 
911         dovertex=kFALSE;
912     }
913     if(dovertex) pvtx=tVertexer.FindPrimaryVertex(esd);
914     if(fDiamondProfile) esd->SetDiamond(fDiamondProfile);
915     
916     if (pvtx)
917     if (pvtx->GetStatus()) {
918        // Store the improved primary vertex
919        esd->SetPrimaryVertex(pvtx);
920        // Propagate the tracks to the DCA to the improved primary vertex
921        Double_t somethingbig = 777.;
922        Double_t bz = esd->GetMagneticField();
923        Int_t nt=esd->GetNumberOfTracks();
924        while (nt--) {
925          AliESDtrack *t = esd->GetTrack(nt);
926          t->RelateToVertex(pvtx, bz, somethingbig);
927        } 
928     }
929
930     if (fRunV0Finder) {
931        // V0 finding
932        AliV0vertexer vtxer;
933        vtxer.Tracks2V0vertices(esd);
934
935        if (fRunCascadeFinder) {
936           // Cascade finding
937           AliCascadeVertexer cvtxer;
938           cvtxer.V0sTracks2CascadeVertices(esd);
939        }
940     }
941  
942     // write ESD
943     if (fCleanESD) CleanESD(esd);
944     if (fWriteESDfriend) {
945       esdf->~AliESDfriend();
946       new (esdf) AliESDfriend(); // Reset...
947       esd->GetESDfriend(esdf);
948     }
949     tree->Fill();
950
951     // write HLT ESD
952     hlttree->Fill();
953
954     if (fCheckPointLevel > 0)  WriteESD(esd, "final"); 
955     esd->Reset();
956     hltesd->Reset();
957     if (fWriteESDfriend) {
958       esdf->~AliESDfriend();
959       new (esdf) AliESDfriend(); // Reset...
960     }
961  
962     gSystem->GetProcInfo(&ProcInfo);
963     AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
964   
965
966   // End of cycle for the in-loop QA
967      if (fRunQA && fInLoopQA) {
968         RunQA(fFillESD.Data(), esd);
969         TString detStr(fFillESD); 
970         for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
971            if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
972            AliQADataMakerRec * qadm = GetQADataMaker(iDet);
973            if (!qadm) continue;
974            qadm->EndOfCycle(AliQA::kRECPOINTS);
975            qadm->EndOfCycle(AliQA::kESDS);
976            qadm->Finish();
977         }
978         if (fRunGlobalQA) {
979            AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
980            if (qadm) {
981               qadm->EndOfCycle(AliQA::kRECPOINTS);
982               qadm->Finish();
983            }
984         }
985      }
986   } 
987   //******** End of the loop over events 
988
989
990
991   tree->GetUserInfo()->Add(esd);
992   hlttree->GetUserInfo()->Add(hltesd);
993   
994   const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();       
995   const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();   
996                  
997    TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());    
998    cdbMapCopy->SetOwner(1);      
999    cdbMapCopy->SetName("cdbMap");        
1000    TIter iter(cdbMap->GetTable());       
1001          
1002    TPair* pair = 0;      
1003    while((pair = dynamic_cast<TPair*> (iter.Next()))){   
1004          TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());   
1005          TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());         
1006          cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));  
1007    }     
1008          
1009    TList *cdbListCopy = new TList();     
1010    cdbListCopy->SetOwner(1);     
1011    cdbListCopy->SetName("cdbList");      
1012          
1013    TIter iter2(cdbList);         
1014          
1015    AliCDBId* id=0;       
1016    while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){         
1017          cdbListCopy->Add(new TObjString(id->ToString().Data()));        
1018    }     
1019          
1020    tree->GetUserInfo()->Add(cdbMapCopy);         
1021    tree->GetUserInfo()->Add(cdbListCopy);
1022
1023
1024   if(fESDPar.Contains("ESD.par")){
1025     AliInfo("Attaching ESD.par to Tree");
1026     TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1027     tree->GetUserInfo()->Add(fn);
1028   }
1029
1030
1031   file->cd();
1032
1033   if (fWriteESDfriend)
1034     tree->SetBranchStatus("ESDfriend*",0);
1035   // we want to have only one tree version number
1036   tree->Write(tree->GetName(),TObject::kOverwrite);
1037   hlttree->Write();
1038
1039   if (fWriteAOD) {
1040     TFile *aodFile = TFile::Open("AliAOD.root", "RECREATE");
1041     ESDFile2AODFile(file, aodFile);
1042     aodFile->Close();
1043   }
1044
1045   gROOT->cd();
1046   CleanUp(file, fileOld);
1047     
1048   // Create tags for the events in the ESD tree (the ESD tree is always present)
1049   // In case of empty events the tags will contain dummy values
1050   AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1051   esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1052   if (fWriteAOD) {
1053     AliAODTagCreator *aodtagCreator = new AliAODTagCreator();
1054     aodtagCreator->CreateAODTags(fFirstEvent,fLastEvent,fGRPList);
1055   }
1056
1057   //Finish QA and end of cycle for out-of-loop QA
1058   if (fRunQA && !fInLoopQA) {
1059      qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS);
1060      //qas.Reset() ;
1061      qas.Run(fRunTracking.Data(), AliQA::kESDS);
1062
1063      if (fRunGlobalQA) {
1064         AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1065         if (qadm) {
1066            qadm->EndOfCycle(AliQA::kRECPOINTS);
1067            qadm->Finish();
1068         }
1069      }
1070   }
1071   
1072   // Cleanup of CDB manager: cache and active storages!
1073   AliCDBManager::Instance()->ClearCache();
1074   
1075   
1076   return kTRUE;
1077 }
1078
1079
1080 //_____________________________________________________________________________
1081 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
1082 {
1083 // run the local reconstruction
1084   static Int_t eventNr=0;
1085   AliCodeTimerAuto("")
1086
1087  //  AliCDBManager* man = AliCDBManager::Instance();
1088 //   Bool_t origCache = man->GetCacheFlag();
1089
1090 //   TString detStr = detectors;
1091 //   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1092 //     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1093 //     AliReconstructor* reconstructor = GetReconstructor(iDet);
1094 //     if (!reconstructor) continue;
1095 //     if (reconstructor->HasLocalReconstruction()) continue;
1096
1097 //     AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1098 //     AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1099     
1100 //     AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));                          
1101 //     AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1102
1103 //     man->SetCacheFlag(kTRUE);
1104 //     TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1105 //     man->GetAll(calibPath); // entries are cached!
1106
1107 //     AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1108      
1109 //     if (fRawReader) {
1110 //       fRawReader->RewindEvents();
1111 //       reconstructor->Reconstruct(fRunLoader, fRawReader);
1112 //     } else {
1113 //       reconstructor->Reconstruct(fRunLoader);
1114 //     }
1115      
1116 //      AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1117     // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1118
1119 //     // unload calibration data
1120 //     man->UnloadFromCache(calibPath);
1121 //     //man->ClearCache();
1122 //   }
1123
1124 //   man->SetCacheFlag(origCache);
1125
1126 //   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1127 //     AliError(Form("the following detectors were not found: %s",
1128 //                   detStr.Data()));
1129 //     if (fStopOnError) return kFALSE;
1130 //   }
1131
1132           eventNr++;
1133   return kTRUE;
1134 }
1135
1136 //_____________________________________________________________________________
1137 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1138 {
1139 // run the local reconstruction
1140
1141   static Int_t eventNr=0;
1142   AliCodeTimerAuto("")
1143
1144   TString detStr = detectors;
1145   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1146     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1147     AliReconstructor* reconstructor = GetReconstructor(iDet);
1148     if (!reconstructor) continue;
1149     AliLoader* loader = fLoader[iDet];
1150     if (!loader) {
1151       AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1152       continue;
1153     }
1154     // conversion of digits
1155     if (fRawReader && reconstructor->HasDigitConversion()) {
1156       AliInfo(Form("converting raw data digits into root objects for %s", 
1157                    fgkDetectorName[iDet]));
1158       AliCodeTimerAuto(Form("converting raw data digits into root objects for %s", 
1159                             fgkDetectorName[iDet]));
1160       loader->LoadDigits("update");
1161       loader->CleanDigits();
1162       loader->MakeDigitsContainer();
1163       TTree* digitsTree = loader->TreeD();
1164       reconstructor->ConvertDigits(fRawReader, digitsTree);
1165       loader->WriteDigits("OVERWRITE");
1166       loader->UnloadDigits();
1167     }
1168     // local reconstruction
1169     AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1170     AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1171     loader->LoadRecPoints("update");
1172     loader->CleanRecPoints();
1173     loader->MakeRecPointsContainer();
1174     TTree* clustersTree = loader->TreeR();
1175     if (fRawReader && !reconstructor->HasDigitConversion()) {
1176       reconstructor->Reconstruct(fRawReader, clustersTree);
1177     } else {
1178       loader->LoadDigits("read");
1179       TTree* digitsTree = loader->TreeD();
1180       if (!digitsTree) {
1181         AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1182         if (fStopOnError) return kFALSE;
1183       } else {
1184         reconstructor->Reconstruct(digitsTree, clustersTree);
1185       }
1186       loader->UnloadDigits();
1187     }
1188
1189     // In-loop QA for local reconstrucion 
1190     if (fRunQA && fInLoopQA) {
1191        AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1192        if (qadm) {
1193           //AliCodeTimerStart
1194           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1195           //AliInfo
1196           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1197
1198           qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1199  
1200           //AliCodeTimerStop
1201           //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1202        }
1203     }
1204
1205     loader->WriteRecPoints("OVERWRITE");
1206     loader->UnloadRecPoints();
1207     AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1208   }
1209
1210   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1211     AliError(Form("the following detectors were not found: %s",
1212                   detStr.Data()));
1213     if (fStopOnError) return kFALSE;
1214   }
1215   eventNr++;
1216   return kTRUE;
1217 }
1218
1219 //_____________________________________________________________________________
1220 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1221 {
1222 // run the barrel tracking
1223
1224   AliCodeTimerAuto("")
1225
1226   AliESDVertex* vertex = NULL;
1227   Double_t vtxPos[3] = {0, 0, 0};
1228   Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1229   TArrayF mcVertex(3); 
1230   if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1231     fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1232     for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1233   }
1234
1235   if (fVertexer) {
1236     if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1237     AliInfo("running the ITS vertex finder");
1238     if (fLoader[0]) fLoader[0]->LoadRecPoints();
1239     vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1240     if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1241     if(!vertex){
1242       AliWarning("Vertex not found");
1243       vertex = new AliESDVertex();
1244       vertex->SetName("default");
1245     }
1246     else {
1247       vertex->SetName("reconstructed");
1248     }
1249
1250   } else {
1251     AliInfo("getting the primary vertex from MC");
1252     vertex = new AliESDVertex(vtxPos, vtxErr);
1253   }
1254
1255   if (vertex) {
1256     vertex->GetXYZ(vtxPos);
1257     vertex->GetSigmaXYZ(vtxErr);
1258   } else {
1259     AliWarning("no vertex reconstructed");
1260     vertex = new AliESDVertex(vtxPos, vtxErr);
1261   }
1262   esd->SetVertex(vertex);
1263   // if SPD multiplicity has been determined, it is stored in the ESD
1264   AliMultiplicity *mult = fVertexer->GetMultiplicity();
1265   if(mult)esd->SetMultiplicity(mult);
1266
1267   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1268     if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1269   }  
1270   delete vertex;
1271
1272   return kTRUE;
1273 }
1274
1275 //_____________________________________________________________________________
1276 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1277 {
1278 // run the HLT barrel tracking
1279
1280   AliCodeTimerAuto("")
1281
1282   if (!fRunLoader) {
1283     AliError("Missing runLoader!");
1284     return kFALSE;
1285   }
1286
1287   AliInfo("running HLT tracking");
1288
1289   // Get a pointer to the HLT reconstructor
1290   AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1291   if (!reconstructor) return kFALSE;
1292
1293   // TPC + ITS
1294   for (Int_t iDet = 1; iDet >= 0; iDet--) {
1295     TString detName = fgkDetectorName[iDet];
1296     AliDebug(1, Form("%s HLT tracking", detName.Data()));
1297     reconstructor->SetOption(detName.Data());
1298     AliTracker *tracker = reconstructor->CreateTracker();
1299     if (!tracker) {
1300       AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1301       if (fStopOnError) return kFALSE;
1302       continue;
1303     }
1304     Double_t vtxPos[3];
1305     Double_t vtxErr[3]={0.005,0.005,0.010};
1306     const AliESDVertex *vertex = esd->GetVertex();
1307     vertex->GetXYZ(vtxPos);
1308     tracker->SetVertex(vtxPos,vtxErr);
1309     if(iDet != 1) {
1310       fLoader[iDet]->LoadRecPoints("read");
1311       TTree* tree = fLoader[iDet]->TreeR();
1312       if (!tree) {
1313         AliError(Form("Can't get the %s cluster tree", detName.Data()));
1314         return kFALSE;
1315       }
1316       tracker->LoadClusters(tree);
1317     }
1318     if (tracker->Clusters2Tracks(esd) != 0) {
1319       AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1320       return kFALSE;
1321     }
1322     if(iDet != 1) {
1323       tracker->UnloadClusters();
1324     }
1325     delete tracker;
1326   }
1327
1328   return kTRUE;
1329 }
1330
1331 //_____________________________________________________________________________
1332 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1333 {
1334 // run the muon spectrometer tracking
1335
1336   AliCodeTimerAuto("")
1337
1338   if (!fRunLoader) {
1339     AliError("Missing runLoader!");
1340     return kFALSE;
1341   }
1342   Int_t iDet = 7; // for MUON
1343
1344   AliInfo("is running...");
1345
1346   // Get a pointer to the MUON reconstructor
1347   AliReconstructor *reconstructor = GetReconstructor(iDet);
1348   if (!reconstructor) return kFALSE;
1349
1350   
1351   TString detName = fgkDetectorName[iDet];
1352   AliDebug(1, Form("%s tracking", detName.Data()));
1353   AliTracker *tracker =  reconstructor->CreateTracker();
1354   if (!tracker) {
1355     AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1356     return kFALSE;
1357   }
1358      
1359   // read RecPoints
1360   fLoader[iDet]->LoadRecPoints("read");  
1361
1362   tracker->LoadClusters(fLoader[iDet]->TreeR());
1363   
1364   Int_t rv = tracker->Clusters2Tracks(esd);
1365   
1366   if ( rv )
1367   {
1368     AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1369     return kFALSE;
1370   }
1371   
1372   fLoader[iDet]->UnloadRecPoints();
1373
1374   tracker->UnloadClusters();
1375   
1376   delete tracker;
1377   
1378   return kTRUE;
1379 }
1380
1381
1382 //_____________________________________________________________________________
1383 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1384 {
1385 // run the barrel tracking
1386   static Int_t eventNr=0;
1387   AliCodeTimerAuto("")
1388
1389   AliInfo("running tracking");
1390
1391   //Fill the ESD with the T0 info (will be used by the TOF) 
1392   if (fReconstructor[11] && fLoader[11]) {
1393     fLoader[11]->LoadRecPoints("READ");
1394     TTree *treeR = fLoader[11]->TreeR();
1395     GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1396   }
1397
1398   // pass 1: TPC + ITS inwards
1399   for (Int_t iDet = 1; iDet >= 0; iDet--) {
1400     if (!fTracker[iDet]) continue;
1401     AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1402
1403     // load clusters
1404     fLoader[iDet]->LoadRecPoints("read");
1405     AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1406     TTree* tree = fLoader[iDet]->TreeR();
1407     if (!tree) {
1408       AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1409       return kFALSE;
1410     }
1411     fTracker[iDet]->LoadClusters(tree);
1412     AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1413     // run tracking
1414     if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1415       AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1416       return kFALSE;
1417     }
1418     if (fCheckPointLevel > 1) {
1419       WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1420     }
1421     // preliminary PID in TPC needed by the ITS tracker
1422     if (iDet == 1) {
1423       GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1424       AliESDpid::MakePID(esd);
1425     } 
1426     AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1427   }
1428
1429   // pass 2: ALL backwards
1430
1431   if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);     
1432
1433   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1434     if (!fTracker[iDet]) continue;
1435     AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1436
1437     // load clusters
1438     if (iDet > 1) {     // all except ITS, TPC
1439       TTree* tree = NULL;
1440       fLoader[iDet]->LoadRecPoints("read");
1441       AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1442       tree = fLoader[iDet]->TreeR();
1443       if (!tree) {
1444         AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1445         return kFALSE;
1446       }
1447       fTracker[iDet]->LoadClusters(tree); 
1448       AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1449     }
1450
1451     // run tracking
1452     if (fTracker[iDet]->PropagateBack(esd) != 0) {
1453       AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1454       //      return kFALSE;
1455     }
1456     if (fCheckPointLevel > 1) {
1457       WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1458     }
1459
1460     // unload clusters
1461     if (iDet > 2) {     // all except ITS, TPC, TRD
1462       fTracker[iDet]->UnloadClusters();
1463       fLoader[iDet]->UnloadRecPoints();
1464     }
1465     // updated PID in TPC needed by the ITS tracker -MI
1466     if (iDet == 1) {
1467       GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1468       AliESDpid::MakePID(esd);
1469     }
1470     AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1471   }
1472
1473   if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);     
1474
1475   // write space-points to the ESD in case alignment data output
1476   // is switched on
1477   if (fWriteAlignmentData)
1478     WriteAlignmentData(esd);
1479
1480   // pass 3: TRD + TPC + ITS refit inwards
1481
1482
1483   for (Int_t iDet = 2; iDet >= 0; iDet--) {
1484     if (!fTracker[iDet]) continue;
1485     AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1486
1487     // run tracking
1488     if (fTracker[iDet]->RefitInward(esd) != 0) {
1489       AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1490       //      return kFALSE;
1491     }
1492     if (fCheckPointLevel > 1) {
1493       WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1494     }
1495     AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1496     // unload clusters
1497     fTracker[iDet]->UnloadClusters();
1498     AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1499     fLoader[iDet]->UnloadRecPoints();
1500     AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1501   }
1502
1503   //
1504   // Propagate track to the vertex - if not done by ITS
1505   //
1506   Int_t ntracks = esd->GetNumberOfTracks();
1507   for (Int_t itrack=0; itrack<ntracks; itrack++){
1508     const Double_t kRadius  = 3;   // beam pipe radius
1509     const Double_t kMaxStep = 5;   // max step
1510     const Double_t kMaxD    = 123456;  // max distance to prim vertex
1511     Double_t       fieldZ   = AliTracker::GetBz();  //
1512     AliESDtrack * track = esd->GetTrack(itrack);
1513     if (!track) continue;
1514     if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1515    AliTracker::PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1516     track->RelateToVertex(esd->GetVertex(),fieldZ, kMaxD);
1517   }
1518   eventNr++;
1519   return kTRUE;
1520 }
1521
1522 //_____________________________________________________________________________
1523 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1524   //
1525   // Remove the data which are not needed for the physics analysis.
1526   //
1527
1528   Int_t nTracks=esd->GetNumberOfTracks();
1529   Int_t nV0s=esd->GetNumberOfV0s();
1530   AliInfo
1531   (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1532
1533   Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1534   Bool_t rc=esd->Clean(cleanPars);
1535
1536   nTracks=esd->GetNumberOfTracks();
1537   nV0s=esd->GetNumberOfV0s();
1538   AliInfo
1539   (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1540
1541   return rc;
1542 }
1543
1544 //_____________________________________________________________________________
1545 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1546 {
1547 // fill the event summary data
1548
1549   AliCodeTimerAuto("")
1550     static Int_t eventNr=0; 
1551   TString detStr = detectors;
1552   
1553   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1554   if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1555     AliReconstructor* reconstructor = GetReconstructor(iDet);
1556     if (!reconstructor) continue;
1557     if (!ReadESD(esd, fgkDetectorName[iDet])) {
1558       AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1559       TTree* clustersTree = NULL;
1560       if (fLoader[iDet]) {
1561         fLoader[iDet]->LoadRecPoints("read");
1562         clustersTree = fLoader[iDet]->TreeR();
1563         if (!clustersTree) {
1564           AliError(Form("Can't get the %s clusters tree", 
1565                         fgkDetectorName[iDet]));
1566           if (fStopOnError) return kFALSE;
1567         }
1568       }
1569       if (fRawReader && !reconstructor->HasDigitConversion()) {
1570         reconstructor->FillESD(fRawReader, clustersTree, esd);
1571       } else {
1572         TTree* digitsTree = NULL;
1573         if (fLoader[iDet]) {
1574           fLoader[iDet]->LoadDigits("read");
1575           digitsTree = fLoader[iDet]->TreeD();
1576           if (!digitsTree) {
1577             AliError(Form("Can't get the %s digits tree", 
1578                           fgkDetectorName[iDet]));
1579             if (fStopOnError) return kFALSE;
1580           }
1581         }
1582         reconstructor->FillESD(digitsTree, clustersTree, esd);
1583         if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1584       }
1585       if (fLoader[iDet]) {
1586         fLoader[iDet]->UnloadRecPoints();
1587       }
1588
1589       if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1590     }
1591   }
1592
1593   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1594     AliError(Form("the following detectors were not found: %s", 
1595                   detStr.Data()));
1596     if (fStopOnError) return kFALSE;
1597   }
1598   AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1599   eventNr++;
1600   return kTRUE;
1601 }
1602
1603 //_____________________________________________________________________________
1604 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1605 {
1606   // Reads the trigger decision which is
1607   // stored in Trigger.root file and fills
1608   // the corresponding esd entries
1609
1610   AliCodeTimerAuto("")
1611   
1612   AliInfo("Filling trigger information into the ESD");
1613
1614   if (fRawReader) {
1615     AliCTPRawStream input(fRawReader);
1616     if (!input.Next()) {
1617       AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1618       return kFALSE;
1619     }
1620     esd->SetTriggerMask(input.GetClassMask());
1621     esd->SetTriggerCluster(input.GetClusterMask());
1622   }
1623   else {
1624     AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1625     if (runloader) {
1626       if (!runloader->LoadTrigger()) {
1627         AliCentralTrigger *aCTP = runloader->GetTrigger();
1628         esd->SetTriggerMask(aCTP->GetClassMask());
1629         esd->SetTriggerCluster(aCTP->GetClusterMask());
1630       }
1631       else {
1632         AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1633         return kFALSE;
1634       }
1635     }
1636     else {
1637       AliError("No run loader is available! The trigger information is not stored in the ESD !");
1638       return kFALSE;
1639     }
1640   }
1641
1642   return kTRUE;
1643 }
1644
1645
1646
1647
1648
1649 //_____________________________________________________________________________
1650 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1651 {
1652   // 
1653   // Filling information from RawReader Header
1654   // 
1655
1656   AliInfo("Filling information from RawReader Header");
1657   esd->SetBunchCrossNumber(0);
1658   esd->SetOrbitNumber(0);
1659   esd->SetPeriodNumber(0);
1660   esd->SetTimeStamp(0);
1661   esd->SetEventType(0);
1662   const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1663   if (eventHeader){
1664
1665     const UInt_t *id = eventHeader->GetP("Id");
1666     esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1667     esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1668     esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1669
1670     esd->SetTimeStamp((eventHeader->Get("Timestamp")));  
1671     esd->SetEventType((eventHeader->Get("Type")));
1672   }
1673
1674   return kTRUE;
1675 }
1676
1677
1678 //_____________________________________________________________________________
1679 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1680 {
1681 // check whether detName is contained in detectors
1682 // if yes, it is removed from detectors
1683
1684   // check if all detectors are selected
1685   if ((detectors.CompareTo("ALL") == 0) ||
1686       detectors.BeginsWith("ALL ") ||
1687       detectors.EndsWith(" ALL") ||
1688       detectors.Contains(" ALL ")) {
1689     detectors = "ALL";
1690     return kTRUE;
1691   }
1692
1693   // search for the given detector
1694   Bool_t result = kFALSE;
1695   if ((detectors.CompareTo(detName) == 0) ||
1696       detectors.BeginsWith(detName+" ") ||
1697       detectors.EndsWith(" "+detName) ||
1698       detectors.Contains(" "+detName+" ")) {
1699     detectors.ReplaceAll(detName, "");
1700     result = kTRUE;
1701   }
1702
1703   // clean up the detectors string
1704   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1705   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1706   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1707
1708   return result;
1709 }
1710
1711 //_____________________________________________________________________________
1712 Bool_t AliReconstruction::InitRunLoader()
1713 {
1714 // get or create the run loader
1715
1716   if (gAlice) delete gAlice;
1717   gAlice = NULL;
1718
1719   if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1720     // load all base libraries to get the loader classes
1721     TString libs = gSystem->GetLibraries();
1722     for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1723       TString detName = fgkDetectorName[iDet];
1724       if (detName == "HLT") continue;
1725       if (libs.Contains("lib" + detName + "base.so")) continue;
1726       gSystem->Load("lib" + detName + "base.so");
1727     }
1728     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1729     if (!fRunLoader) {
1730       AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1731       CleanUp();
1732       return kFALSE;
1733     }
1734     fRunLoader->CdGAFile();
1735     if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1736       if (fRunLoader->LoadgAlice() == 0) {
1737         gAlice = fRunLoader->GetAliRun();
1738         AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1739       }
1740     }
1741     if (!gAlice && !fRawReader) {
1742       AliError(Form("no gAlice object found in file %s",
1743                     fGAliceFileName.Data()));
1744       CleanUp();
1745       return kFALSE;
1746     }
1747
1748     //PH This is a temporary fix to give access to the kinematics
1749     //PH that is needed for the labels of ITS clusters
1750     fRunLoader->LoadHeader();
1751     fRunLoader->LoadKinematics();
1752
1753   } else {               // galice.root does not exist
1754     if (!fRawReader) {
1755       AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1756       CleanUp();
1757       return kFALSE;
1758     }
1759     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1760                                     AliConfig::GetDefaultEventFolderName(),
1761                                     "recreate");
1762     if (!fRunLoader) {
1763       AliError(Form("could not create run loader in file %s", 
1764                     fGAliceFileName.Data()));
1765       CleanUp();
1766       return kFALSE;
1767     }
1768     fRunLoader->MakeTree("E");
1769     Int_t iEvent = 0;
1770     while (fRawReader->NextEvent()) {
1771       fRunLoader->SetEventNumber(iEvent);
1772       fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(), 
1773                                      iEvent, iEvent);
1774       fRunLoader->MakeTree("H");
1775       fRunLoader->TreeE()->Fill();
1776       iEvent++;
1777     }
1778     fRawReader->RewindEvents();
1779     if (fNumberOfEventsPerFile > 0)
1780       fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1781     else
1782       fRunLoader->SetNumberOfEventsPerFile(iEvent);
1783     fRunLoader->WriteHeader("OVERWRITE");
1784     fRunLoader->CdGAFile();
1785     fRunLoader->Write(0, TObject::kOverwrite);
1786 //    AliTracker::SetFieldMap(???);
1787   }
1788
1789   return kTRUE;
1790 }
1791
1792 //_____________________________________________________________________________
1793 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1794 {
1795 // get the reconstructor object and the loader for a detector
1796
1797   if (fReconstructor[iDet]) return fReconstructor[iDet];
1798
1799   // load the reconstructor object
1800   TPluginManager* pluginManager = gROOT->GetPluginManager();
1801   TString detName = fgkDetectorName[iDet];
1802   TString recName = "Ali" + detName + "Reconstructor";
1803   if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1804
1805   AliReconstructor* reconstructor = NULL;
1806   // first check if a plugin is defined for the reconstructor
1807   TPluginHandler* pluginHandler = 
1808     pluginManager->FindHandler("AliReconstructor", detName);
1809   // if not, add a plugin for it
1810   if (!pluginHandler) {
1811     AliDebug(1, Form("defining plugin for %s", recName.Data()));
1812     TString libs = gSystem->GetLibraries();
1813     if (libs.Contains("lib" + detName + "base.so") ||
1814         (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1815       pluginManager->AddHandler("AliReconstructor", detName, 
1816                                 recName, detName + "rec", recName + "()");
1817     } else {
1818       pluginManager->AddHandler("AliReconstructor", detName, 
1819                                 recName, detName, recName + "()");
1820     }
1821     pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1822   }
1823   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1824     reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1825   }
1826   if (reconstructor) {
1827     TObject* obj = fOptions.FindObject(detName.Data());
1828     if (obj) reconstructor->SetOption(obj->GetTitle());
1829     reconstructor->Init();
1830     fReconstructor[iDet] = reconstructor;
1831   }
1832
1833   // get or create the loader
1834   if (detName != "HLT") {
1835     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1836     if (!fLoader[iDet]) {
1837       AliConfig::Instance()
1838         ->CreateDetectorFolders(fRunLoader->GetEventFolder(), 
1839                                 detName, detName);
1840       // first check if a plugin is defined for the loader
1841       pluginHandler = 
1842         pluginManager->FindHandler("AliLoader", detName);
1843       // if not, add a plugin for it
1844       if (!pluginHandler) {
1845         TString loaderName = "Ali" + detName + "Loader";
1846         AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1847         pluginManager->AddHandler("AliLoader", detName, 
1848                                   loaderName, detName + "base", 
1849                                   loaderName + "(const char*, TFolder*)");
1850         pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1851       }
1852       if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1853         fLoader[iDet] = 
1854           (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(), 
1855                                                  fRunLoader->GetEventFolder());
1856       }
1857       if (!fLoader[iDet]) {   // use default loader
1858         fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1859       }
1860       if (!fLoader[iDet]) {
1861         AliWarning(Form("couldn't get loader for %s", detName.Data()));
1862         if (fStopOnError) return NULL;
1863       } else {
1864         fRunLoader->AddLoader(fLoader[iDet]);
1865         fRunLoader->CdGAFile();
1866         if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1867         fRunLoader->Write(0, TObject::kOverwrite);
1868       }
1869     }
1870   }
1871       
1872   return reconstructor;
1873 }
1874
1875 //_____________________________________________________________________________
1876 Bool_t AliReconstruction::CreateVertexer()
1877 {
1878 // create the vertexer
1879
1880   fVertexer = NULL;
1881   AliReconstructor* itsReconstructor = GetReconstructor(0);
1882   if (itsReconstructor) {
1883     fVertexer = itsReconstructor->CreateVertexer();
1884   }
1885   if (!fVertexer) {
1886     AliWarning("couldn't create a vertexer for ITS");
1887     if (fStopOnError) return kFALSE;
1888   }
1889
1890   return kTRUE;
1891 }
1892
1893 //_____________________________________________________________________________
1894 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
1895 {
1896 // create the trackers
1897
1898   TString detStr = detectors;
1899   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1900     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1901     AliReconstructor* reconstructor = GetReconstructor(iDet);
1902     if (!reconstructor) continue;
1903     TString detName = fgkDetectorName[iDet];
1904     if (detName == "HLT") {
1905       fRunHLTTracking = kTRUE;
1906       continue;
1907     }
1908     if (detName == "MUON") {
1909       fRunMuonTracking = kTRUE;
1910       continue;
1911     }
1912
1913
1914     fTracker[iDet] = reconstructor->CreateTracker();
1915     if (!fTracker[iDet] && (iDet < 7)) {
1916       AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1917       if (fStopOnError) return kFALSE;
1918     }
1919     AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
1920   }
1921
1922   return kTRUE;
1923 }
1924
1925 //_____________________________________________________________________________
1926 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
1927 {
1928 // delete trackers and the run loader and close and delete the file
1929
1930   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1931     delete fReconstructor[iDet];
1932     fReconstructor[iDet] = NULL;
1933     fLoader[iDet] = NULL;
1934     delete fTracker[iDet];
1935     fTracker[iDet] = NULL;
1936 //    delete fQADataMaker[iDet];
1937 //    fQADataMaker[iDet] = NULL;
1938   }
1939   delete fVertexer;
1940   fVertexer = NULL;
1941   
1942   if(!(AliCDBManager::Instance()->GetCacheFlag())) {
1943         delete fDiamondProfile;
1944         fDiamondProfile = NULL;
1945   }
1946
1947   delete fGRPList;
1948   fGRPList = NULL;
1949
1950   delete fRunLoader;
1951   fRunLoader = NULL;
1952   delete fRawReader;
1953   fRawReader = NULL;
1954
1955   if (file) {
1956     file->Close();
1957     delete file;
1958   }
1959
1960   if (fileOld) {
1961     fileOld->Close();
1962     delete fileOld;
1963     gSystem->Unlink("AliESDs.old.root");
1964   }
1965 }
1966
1967 //_____________________________________________________________________________
1968
1969 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
1970 {
1971 // read the ESD event from a file
1972
1973   if (!esd) return kFALSE;
1974   char fileName[256];
1975   sprintf(fileName, "ESD_%d.%d_%s.root", 
1976           esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
1977   if (gSystem->AccessPathName(fileName)) return kFALSE;
1978
1979   AliInfo(Form("reading ESD from file %s", fileName));
1980   AliDebug(1, Form("reading ESD from file %s", fileName));
1981   TFile* file = TFile::Open(fileName);
1982   if (!file || !file->IsOpen()) {
1983     AliError(Form("opening %s failed", fileName));
1984     delete file;
1985     return kFALSE;
1986   }
1987
1988   gROOT->cd();
1989   delete esd;
1990   esd = (AliESDEvent*) file->Get("ESD");
1991   file->Close();
1992   delete file;
1993   return kTRUE;
1994
1995 }
1996
1997
1998
1999 //_____________________________________________________________________________
2000 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2001 {
2002 // write the ESD event to a file
2003
2004   if (!esd) return;
2005   char fileName[256];
2006   sprintf(fileName, "ESD_%d.%d_%s.root", 
2007           esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2008
2009   AliDebug(1, Form("writing ESD to file %s", fileName));
2010   TFile* file = TFile::Open(fileName, "recreate");
2011   if (!file || !file->IsOpen()) {
2012     AliError(Form("opening %s failed", fileName));
2013   } else {
2014     esd->Write("ESD");
2015     file->Close();
2016   }
2017   delete file;
2018 }
2019
2020
2021
2022
2023
2024 //_____________________________________________________________________________
2025 void AliReconstruction::ESDFile2AODFile(TFile* esdFile, TFile* aodFile)
2026 {
2027   // write all files from the given esd file to an aod file
2028
2029   // create an AliAOD object 
2030   AliAODEvent *aod = new AliAODEvent();
2031   aod->CreateStdContent();
2032   
2033   // go to the file
2034   aodFile->cd();
2035   
2036   // create the tree
2037   TTree *aodTree = new TTree("aodTree", "AliAOD tree");
2038   aodTree->Branch(aod->GetList());
2039
2040   // connect to ESD
2041   TTree *t = (TTree*) esdFile->Get("esdTree");
2042   AliESDEvent *esd = new AliESDEvent();
2043   esd->ReadFromTree(t);
2044
2045   Int_t nEvents = t->GetEntries();
2046
2047   // set arrays and pointers
2048   Float_t posF[3];
2049   Double_t pos[3];
2050   Double_t p[3];
2051   Double_t p_pos[3];
2052   Double_t p_neg[3];
2053   Double_t covVtx[6];
2054   Double_t covTr[21];
2055   Double_t pid[10];
2056
2057   // loop over events and fill them
2058   for (Int_t iEvent = 0; iEvent < nEvents; ++iEvent) {
2059     //cout << "event: " << iEvent << endl;
2060     t->GetEntry(iEvent);
2061
2062     // Multiplicity information needed by the header (to be revised!)
2063     Int_t nTracks   = esd->GetNumberOfTracks();
2064     Int_t nPosTracks = 0;
2065     for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) 
2066       if (esd->GetTrack(iTrack)->Charge()> 0) nPosTracks++;
2067
2068     // Access the header
2069     AliAODHeader *header = aod->GetHeader();
2070
2071     // fill the header
2072     header->SetRunNumber       (esd->GetRunNumber()       );
2073     header->SetBunchCrossNumber(esd->GetBunchCrossNumber());
2074     header->SetOrbitNumber     (esd->GetOrbitNumber()     );
2075     header->SetPeriodNumber    (esd->GetPeriodNumber()    );
2076     header->SetTriggerMask     (esd->GetTriggerMask()     ); 
2077     header->SetTriggerCluster  (esd->GetTriggerCluster()  );
2078     header->SetEventType       (esd->GetEventType()       );
2079     header->SetMagneticField   (esd->GetMagneticField()   );
2080     header->SetZDCN1Energy     (esd->GetZDCN1Energy()     );
2081     header->SetZDCP1Energy     (esd->GetZDCP1Energy()     );
2082     header->SetZDCN2Energy     (esd->GetZDCN2Energy()     );
2083     header->SetZDCP2Energy     (esd->GetZDCP2Energy()     );
2084     header->SetZDCEMEnergy     (esd->GetZDCEMEnergy(0),esd->GetZDCEMEnergy(1));
2085     header->SetRefMultiplicity   (nTracks);
2086     header->SetRefMultiplicityPos(nPosTracks);
2087     header->SetRefMultiplicityNeg(nTracks - nPosTracks);
2088     header->SetMuonMagFieldScale(-999.); // FIXME
2089     header->SetCentrality(-999.);        // FIXME
2090
2091     Int_t nV0s      = esd->GetNumberOfV0s();
2092     Int_t nCascades = esd->GetNumberOfCascades();
2093     Int_t nKinks    = esd->GetNumberOfKinks();
2094         Int_t nVertices = nV0s + 2*nCascades /*could lead to two vertices, one V0 and the Xi */+ nKinks + 1 /* = prim. vtx*/;    
2095         Int_t nJets     = 0;
2096     Int_t nCaloClus = esd->GetNumberOfCaloClusters();
2097     Int_t nFmdClus  = 0;
2098     Int_t nPmdClus  = esd->GetNumberOfPmdTracks();
2099    
2100     aod->ResetStd(nTracks, nVertices, nV0s+nCascades, nJets, nCaloClus, nFmdClus, nPmdClus);
2101     
2102     // Array to take into account the tracks already added to the AOD
2103     Bool_t * usedTrack = NULL;
2104     if (nTracks>0) {
2105       usedTrack = new Bool_t[nTracks];
2106       for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) usedTrack[iTrack]=kFALSE;
2107     }
2108     // Array to take into account the V0s already added to the AOD
2109     Bool_t * usedV0 = NULL;
2110     if (nV0s>0) {
2111       usedV0 = new Bool_t[nV0s];
2112       for (Int_t iV0=0; iV0<nV0s; ++iV0) usedV0[iV0]=kFALSE;
2113     }
2114     // Array to take into account the kinks already added to the AOD
2115     Bool_t * usedKink = NULL;
2116     if (nKinks>0) {
2117       usedKink = new Bool_t[nKinks];
2118       for (Int_t iKink=0; iKink<nKinks; ++iKink) usedKink[iKink]=kFALSE;
2119     }
2120     
2121     // Access to the AOD container of vertices
2122     TClonesArray &vertices = *(aod->GetVertices());
2123     Int_t jVertices=0;
2124
2125     // Access to the AOD container of tracks
2126     TClonesArray &tracks = *(aod->GetTracks());
2127     Int_t jTracks=0; 
2128    
2129     // Access to the AOD container of V0s
2130     TClonesArray &V0s = *(aod->GetV0s());
2131     Int_t jV0s=0;
2132     
2133     // Add primary vertex. The primary tracks will be defined
2134     // after the loops on the composite objects (V0, cascades, kinks)
2135     const AliESDVertex *vtx = esd->GetPrimaryVertex();
2136       
2137     vtx->GetXYZ(pos); // position
2138     vtx->GetCovMatrix(covVtx); //covariance matrix
2139
2140     AliAODVertex * primary = new(vertices[jVertices++])
2141       AliAODVertex(pos, covVtx, vtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
2142          
2143
2144     AliAODTrack *aodTrack = 0x0;
2145     
2146     // Create vertices starting from the most complex objects
2147
2148     // Cascades
2149     for (Int_t nCascade = 0; nCascade < nCascades; ++nCascade) {
2150       AliESDcascade *cascade = esd->GetCascade(nCascade);
2151       
2152       cascade->GetXYZ(pos[0], pos[1], pos[2]);
2153       cascade->GetPosCovXi(covVtx);
2154      
2155       // Add the cascade vertex
2156       AliAODVertex * vcascade = new(vertices[jVertices++]) AliAODVertex(pos,
2157                                                                         covVtx,
2158                                                                         cascade->GetChi2Xi(), // = chi2/NDF since NDF = 2*2-3
2159                                                                         primary,
2160                                                                         nCascade,
2161                                                                         AliAODVertex::kCascade);
2162
2163       primary->AddDaughter(vcascade); // the cascade 'particle' (represented by a vertex) is added as a daughter to the primary vertex
2164
2165       // Add the V0 from the cascade. The ESD class have to be optimized...
2166       // Now we have to search for the corresponding V0 in the list of V0s
2167       // using the indeces of the positive and negative tracks
2168
2169       Int_t posFromV0 = cascade->GetPindex();
2170       Int_t negFromV0 = cascade->GetNindex();
2171
2172
2173       AliESDv0 * v0 = 0x0;
2174       Int_t indV0 = -1;
2175
2176       for (Int_t iV0=0; iV0<nV0s; ++iV0) {
2177
2178         v0 = esd->GetV0(iV0);
2179         Int_t posV0 = v0->GetPindex();
2180         Int_t negV0 = v0->GetNindex();
2181
2182         if (posV0==posFromV0 && negV0==negFromV0) {
2183           indV0 = iV0;
2184           break;
2185         }
2186       }
2187
2188       AliAODVertex * vV0FromCascade = 0x0;
2189
2190       if (indV0>-1 && !usedV0[indV0]) {
2191         
2192         // the V0 exists in the array of V0s and is not used
2193
2194         usedV0[indV0] = kTRUE;
2195         
2196         v0->GetXYZ(pos[0], pos[1], pos[2]);
2197         v0->GetPosCov(covVtx);
2198         
2199         vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2200                                                                  covVtx,
2201                                                                  v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2202                                                                  vcascade,
2203                                                                  indV0,
2204                                                                  AliAODVertex::kV0);
2205       } else {
2206
2207         // the V0 doesn't exist in the array of V0s or was used
2208         cerr << "Error: event " << iEvent << " cascade " << nCascade
2209              << " The V0 " << indV0 
2210              << " doesn't exist in the array of V0s or was used!" << endl;
2211
2212         cascade->GetXYZ(pos[0], pos[1], pos[2]);
2213         cascade->GetPosCov(covVtx);
2214       
2215         vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2216                                                                  covVtx,
2217                                                                  v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2218                                                                  vcascade,
2219                                                                  indV0,
2220                                                                  AliAODVertex::kV0);
2221         vcascade->AddDaughter(vV0FromCascade);
2222
2223       }
2224
2225       // Add the positive tracks from the V0
2226
2227       if (! usedTrack[posFromV0]) {
2228
2229         usedTrack[posFromV0] = kTRUE;
2230
2231         AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2232         esdTrack->GetPxPyPz(p_pos);
2233         esdTrack->GetXYZ(pos);
2234         esdTrack->GetCovarianceXYZPxPyPz(covTr);
2235         esdTrack->GetESDpid(pid);
2236         
2237         vV0FromCascade->AddDaughter(aodTrack =
2238                                     new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2239                                            esdTrack->GetLabel(), 
2240                                            p_pos, 
2241                                            kTRUE,
2242                                            pos,
2243                                            kFALSE,
2244                                            covTr, 
2245                                            (Short_t)esdTrack->Charge(),
2246                                            esdTrack->GetITSClusterMap(), 
2247                                            pid,
2248                                            vV0FromCascade,
2249                                            kTRUE,  // check if this is right
2250                                            kFALSE, // check if this is right
2251                                            AliAODTrack::kSecondary)
2252                 );
2253         aodTrack->ConvertAliPIDtoAODPID();
2254       }
2255       else {
2256         cerr << "Error: event " << iEvent << " cascade " << nCascade
2257              << " track " << posFromV0 << " has already been used!" << endl;
2258       }
2259
2260       // Add the negative tracks from the V0
2261
2262       if (!usedTrack[negFromV0]) {
2263         
2264         usedTrack[negFromV0] = kTRUE;
2265         
2266         AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2267         esdTrack->GetPxPyPz(p_neg);
2268         esdTrack->GetXYZ(pos);
2269         esdTrack->GetCovarianceXYZPxPyPz(covTr);
2270         esdTrack->GetESDpid(pid);
2271         
2272         vV0FromCascade->AddDaughter(aodTrack =
2273                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2274                                            esdTrack->GetLabel(),
2275                                            p_neg,
2276                                            kTRUE,
2277                                            pos,
2278                                            kFALSE,
2279                                            covTr, 
2280                                            (Short_t)esdTrack->Charge(),
2281                                            esdTrack->GetITSClusterMap(), 
2282                                            pid,
2283                                            vV0FromCascade,
2284                                            kTRUE,  // check if this is right
2285                                            kFALSE, // check if this is right
2286                                            AliAODTrack::kSecondary)
2287                 );
2288         aodTrack->ConvertAliPIDtoAODPID();
2289       }
2290       else {
2291         cerr << "Error: event " << iEvent << " cascade " << nCascade
2292              << " track " << negFromV0 << " has already been used!" << endl;
2293       }
2294
2295       // add it to the V0 array as well
2296       Double_t d0[2] = { -999., -99.};
2297       // counting is probably wrong
2298       new(V0s[jV0s++]) AliAODv0(vV0FromCascade, -999., -99., p_pos, p_neg, d0); // to be refined
2299
2300       // Add the bachelor track from the cascade
2301
2302       Int_t bachelor = cascade->GetBindex();
2303       
2304       if(!usedTrack[bachelor]) {
2305       
2306         usedTrack[bachelor] = kTRUE;
2307         
2308         AliESDtrack *esdTrack = esd->GetTrack(bachelor);
2309         esdTrack->GetPxPyPz(p);
2310         esdTrack->GetXYZ(pos);
2311         esdTrack->GetCovarianceXYZPxPyPz(covTr);
2312         esdTrack->GetESDpid(pid);
2313
2314         vcascade->AddDaughter(aodTrack =
2315                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2316                                            esdTrack->GetLabel(),
2317                                            p,
2318                                            kTRUE,
2319                                            pos,
2320                                            kFALSE,
2321                                            covTr, 
2322                                            (Short_t)esdTrack->Charge(),
2323                                            esdTrack->GetITSClusterMap(), 
2324                                            pid,
2325                                            vcascade,
2326                                            kTRUE,  // check if this is right
2327                                            kFALSE, // check if this is right
2328                                            AliAODTrack::kSecondary)
2329                 );
2330         aodTrack->ConvertAliPIDtoAODPID();
2331      }
2332       else {
2333         cerr << "Error: event " << iEvent << " cascade " << nCascade
2334              << " track " << bachelor << " has already been used!" << endl;
2335       }
2336       
2337       // Add the primary track of the cascade (if any)
2338       
2339     } // end of the loop on cascades
2340  
2341     // V0s
2342         
2343     for (Int_t nV0 = 0; nV0 < nV0s; ++nV0) {
2344
2345       if (usedV0[nV0]) continue; // skip if aready added to the AOD
2346
2347       AliESDv0 *v0 = esd->GetV0(nV0); 
2348      
2349       v0->GetXYZ(pos[0], pos[1], pos[2]);
2350       v0->GetPosCov(covVtx);
2351
2352       AliAODVertex * vV0 = 
2353         new(vertices[jVertices++]) AliAODVertex(pos,
2354                                                 covVtx,
2355                                                 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2356                                                 primary,
2357                                                 nV0,
2358                                                 AliAODVertex::kV0);
2359       primary->AddDaughter(vV0);
2360
2361       Int_t posFromV0 = v0->GetPindex();
2362       Int_t negFromV0 = v0->GetNindex();
2363       
2364       // Add the positive tracks from the V0
2365
2366       if (!usedTrack[posFromV0]) {
2367         
2368         usedTrack[posFromV0] = kTRUE;
2369
2370         AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2371         esdTrack->GetPxPyPz(p_pos);
2372         esdTrack->GetXYZ(pos);
2373         esdTrack->GetCovarianceXYZPxPyPz(covTr);
2374         esdTrack->GetESDpid(pid);
2375         
2376         vV0->AddDaughter(aodTrack =
2377                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2378                                            esdTrack->GetLabel(), 
2379                                            p_pos, 
2380                                            kTRUE,
2381                                            pos,
2382                                            kFALSE,
2383                                            covTr, 
2384                                            (Short_t)esdTrack->Charge(),
2385                                            esdTrack->GetITSClusterMap(), 
2386                                            pid,
2387                                            vV0,
2388                                            kTRUE,  // check if this is right
2389                                            kFALSE, // check if this is right
2390                                            AliAODTrack::kSecondary)
2391                 );
2392         aodTrack->ConvertAliPIDtoAODPID();
2393       }
2394       else {
2395         cerr << "Error: event " << iEvent << " V0 " << nV0
2396              << " track " << posFromV0 << " has already been used!" << endl;
2397       }
2398
2399       // Add the negative tracks from the V0
2400
2401       if (!usedTrack[negFromV0]) {
2402
2403         usedTrack[negFromV0] = kTRUE;
2404
2405         AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2406         esdTrack->GetPxPyPz(p_neg);
2407         esdTrack->GetXYZ(pos);
2408         esdTrack->GetCovarianceXYZPxPyPz(covTr);
2409         esdTrack->GetESDpid(pid);
2410
2411         vV0->AddDaughter(aodTrack =
2412                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2413                                            esdTrack->GetLabel(),
2414                                            p_neg,
2415                                            kTRUE,
2416                                            pos,
2417                                            kFALSE,
2418                                            covTr, 
2419                                            (Short_t)esdTrack->Charge(),
2420                                            esdTrack->GetITSClusterMap(), 
2421                                            pid,
2422                                            vV0,
2423                                            kTRUE,  // check if this is right
2424                                            kFALSE, // check if this is right
2425                                            AliAODTrack::kSecondary)
2426                 );
2427         aodTrack->ConvertAliPIDtoAODPID();
2428       }
2429       else {
2430         cerr << "Error: event " << iEvent << " V0 " << nV0
2431              << " track " << negFromV0 << " has already been used!" << endl;
2432       }
2433
2434       // add it to the V0 array as well
2435       Double_t d0[2] = { 999., 99.};
2436       new(V0s[jV0s++]) AliAODv0(vV0, 999., 99., p_pos, p_neg, d0); // to be refined
2437     }
2438         V0s.Expand(jV0s);        
2439     // end of the loop on V0s
2440     
2441     // Kinks: it is a big mess the access to the information in the kinks
2442     // The loop is on the tracks in order to find the mother and daugther of each kink
2443
2444
2445     for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) {
2446
2447       AliESDtrack * esdTrack = esd->GetTrack(iTrack);
2448
2449       Int_t ikink = esdTrack->GetKinkIndex(0);
2450
2451       if (ikink) {
2452         // Negative kink index: mother, positive: daughter
2453
2454         // Search for the second track of the kink
2455
2456         for (Int_t jTrack = iTrack+1; jTrack<nTracks; ++jTrack) {
2457
2458           AliESDtrack * esdTrack1 = esd->GetTrack(jTrack);
2459
2460           Int_t jkink = esdTrack1->GetKinkIndex(0);
2461
2462           if ( TMath::Abs(ikink)==TMath::Abs(jkink) ) {
2463
2464             // The two tracks are from the same kink
2465           
2466             if (usedKink[TMath::Abs(ikink)-1]) continue; // skip used kinks
2467
2468             Int_t imother = -1;
2469             Int_t idaughter = -1;
2470
2471             if (ikink<0 && jkink>0) {
2472
2473               imother = iTrack;
2474               idaughter = jTrack;
2475             }
2476             else if (ikink>0 && jkink<0) {
2477
2478               imother = jTrack;
2479               idaughter = iTrack;
2480             }
2481             else {
2482               cerr << "Error: Wrong combination of kink indexes: "
2483               << ikink << " " << jkink << endl;
2484               continue;
2485             }
2486
2487             // Add the mother track
2488
2489             AliAODTrack * mother = NULL;
2490
2491             if (!usedTrack[imother]) {
2492         
2493               usedTrack[imother] = kTRUE;
2494         
2495               AliESDtrack *esdTrack = esd->GetTrack(imother);
2496               esdTrack->GetPxPyPz(p);
2497               esdTrack->GetXYZ(pos);
2498               esdTrack->GetCovarianceXYZPxPyPz(covTr);
2499               esdTrack->GetESDpid(pid);
2500
2501               mother = 
2502                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2503                                            esdTrack->GetLabel(),
2504                                            p,
2505                                            kTRUE,
2506                                            pos,
2507                                            kFALSE,
2508                                            covTr, 
2509                                            (Short_t)esdTrack->Charge(),
2510                                            esdTrack->GetITSClusterMap(), 
2511                                            pid,
2512                                            primary,
2513                                            kTRUE, // check if this is right
2514                                            kTRUE, // check if this is right
2515                                            AliAODTrack::kPrimary);
2516               primary->AddDaughter(mother);
2517               mother->ConvertAliPIDtoAODPID();
2518             }
2519             else {
2520               cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2521               << " track " << imother << " has already been used!" << endl;
2522             }
2523
2524             // Add the kink vertex
2525             AliESDkink * kink = esd->GetKink(TMath::Abs(ikink)-1);
2526
2527             AliAODVertex * vkink = 
2528             new(vertices[jVertices++]) AliAODVertex(kink->GetPosition(),
2529                                                     NULL,
2530                                                     0.,
2531                                                     mother,
2532                                                     esdTrack->GetID(), // This is the track ID of the mother's track!
2533                                                     AliAODVertex::kKink);
2534             // Add the daughter track
2535
2536             AliAODTrack * daughter = NULL;
2537
2538             if (!usedTrack[idaughter]) {
2539         
2540               usedTrack[idaughter] = kTRUE;
2541         
2542               AliESDtrack *esdTrack = esd->GetTrack(idaughter);
2543               esdTrack->GetPxPyPz(p);
2544               esdTrack->GetXYZ(pos);
2545               esdTrack->GetCovarianceXYZPxPyPz(covTr);
2546               esdTrack->GetESDpid(pid);
2547
2548               daughter = 
2549                 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2550                                            esdTrack->GetLabel(),
2551                                            p,
2552                                            kTRUE,
2553                                            pos,
2554                                            kFALSE,
2555                                            covTr, 
2556                                            (Short_t)esdTrack->Charge(),
2557                                            esdTrack->GetITSClusterMap(), 
2558                                            pid,
2559                                            vkink,
2560                                            kTRUE, // check if this is right
2561                                            kTRUE, // check if this is right
2562                                            AliAODTrack::kPrimary);
2563               vkink->AddDaughter(daughter);
2564               daughter->ConvertAliPIDtoAODPID();
2565             }
2566             else {
2567               cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2568               << " track " << idaughter << " has already been used!" << endl;
2569             }
2570           }
2571         }
2572       }
2573     }
2574     vertices.Expand(jVertices);
2575
2576     // Tracks (primary and orphan)
2577     for (Int_t nTrack = 0; nTrack < nTracks; ++nTrack) {
2578
2579       if (usedTrack[nTrack]) continue;
2580
2581       AliESDtrack *esdTrack = esd->GetTrack(nTrack);
2582       esdTrack->GetPxPyPz(p);
2583       esdTrack->GetXYZ(pos);
2584       esdTrack->GetCovarianceXYZPxPyPz(covTr);
2585       esdTrack->GetESDpid(pid);
2586
2587       Float_t impactXY, impactZ;
2588
2589       esdTrack->GetImpactParameters(impactXY,impactZ);
2590
2591       if (impactXY<3.) {
2592         // track inside the beam pipe
2593       
2594         primary->AddDaughter(aodTrack =
2595             new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2596                                          esdTrack->GetLabel(),
2597                                          p,
2598                                          kTRUE,
2599                                          pos,
2600                                          kFALSE,
2601                                          covTr, 
2602                                          (Short_t)esdTrack->Charge(),
2603                                          esdTrack->GetITSClusterMap(), 
2604                                          pid,
2605                                          primary,
2606                                          kTRUE, // check if this is right
2607                                          kTRUE, // check if this is right
2608                                          AliAODTrack::kPrimary)
2609             );
2610         aodTrack->ConvertAliPIDtoAODPID();
2611       }
2612       else {
2613         // outside the beam pipe: orphan track
2614         // Don't write them anymore!
2615         continue;
2616       } 
2617     } // end of loop on tracks
2618     
2619     // muon tracks
2620     Int_t nMuTracks = esd->GetNumberOfMuonTracks();
2621     for (Int_t nMuTrack = 0; nMuTrack < nMuTracks; ++nMuTrack) {
2622       
2623       AliESDMuonTrack *esdMuTrack = esd->GetMuonTrack(nMuTrack);     
2624       p[0] = esdMuTrack->Px(); 
2625       p[1] = esdMuTrack->Py(); 
2626       p[2] = esdMuTrack->Pz();
2627       pos[0] = primary->GetX(); 
2628       pos[1] = primary->GetY(); 
2629       pos[2] = primary->GetZ();
2630       
2631       // has to be changed once the muon pid is provided by the ESD
2632       for (Int_t i = 0; i < 10; pid[i++] = 0.); pid[AliAODTrack::kMuon]=1.;
2633       
2634       primary->AddDaughter(aodTrack =
2635           new(tracks[jTracks++]) AliAODTrack(0, // no ID provided
2636                                              0, // no label provided
2637                                              p,
2638                                              kTRUE,
2639                                              pos,
2640                                              kFALSE,
2641                                              NULL, // no covariance matrix provided
2642                                              esdMuTrack->Charge(),
2643                                              0, // ITSClusterMap is set below
2644                                              pid,
2645                                              primary,
2646                                              kFALSE,  // muon tracks are not used to fit the primary vtx
2647                                              kFALSE,  // not used for vertex fit
2648                                              AliAODTrack::kPrimary)
2649           );
2650
2651       aodTrack->SetHitsPatternInTrigCh(esdMuTrack->GetHitsPatternInTrigCh());
2652       Int_t track2Trigger = esdMuTrack->GetMatchTrigger();
2653       aodTrack->SetMatchTrigger(track2Trigger);
2654       if (track2Trigger) 
2655         aodTrack->SetChi2MatchTrigger(esdMuTrack->GetChi2MatchTrigger());
2656       else 
2657         aodTrack->SetChi2MatchTrigger(0.);
2658     }
2659     tracks.Expand(jTracks); // remove 'empty slots' due to unwritten tracks
2660         
2661     // Access to the AOD container of PMD clusters
2662     TClonesArray &pmdClusters = *(aod->GetPmdClusters());
2663     Int_t jPmdClusters=0;
2664   
2665     for (Int_t iPmd = 0; iPmd < nPmdClus; ++iPmd) {
2666       // file pmd clusters, to be revised!
2667       AliESDPmdTrack *pmdTrack = esd->GetPmdTrack(iPmd);
2668       Int_t nLabel = 0;
2669       Int_t *label = 0x0;
2670       Double_t pos[3] = { pmdTrack->GetClusterX(), pmdTrack->GetClusterY(), pmdTrack->GetClusterZ() };
2671       Double_t pid[9] = { 0., 0., 0., 0., 0., 0., 0., 0., 0. }; // to be revised!
2672       // type not set!
2673       // assoc cluster not set
2674       new(pmdClusters[jPmdClusters++]) AliAODPmdCluster(iPmd, nLabel, label, pmdTrack->GetClusterADC(), pos, pid);
2675     }
2676
2677     // Access to the AOD container of clusters
2678     TClonesArray &caloClusters = *(aod->GetCaloClusters());
2679     Int_t jClusters=0;
2680  
2681     for (Int_t iClust=0; iClust<nCaloClus; ++iClust) {
2682
2683       AliESDCaloCluster * cluster = esd->GetCaloCluster(iClust);
2684
2685       Int_t id = cluster->GetID();
2686       Int_t nLabel = 0;
2687       Int_t *label = 0x0;
2688       Float_t energy = cluster->E();
2689       cluster->GetPosition(posF);
2690       Char_t ttype=AliAODCluster::kUndef;
2691
2692       if (cluster->GetClusterType() == AliESDCaloCluster::kPHOSCluster) {
2693         ttype=AliAODCluster::kPHOSNeutral;
2694       } 
2695       else if (cluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {
2696         ttype = AliAODCluster::kEMCALClusterv1;
2697       }
2698
2699       
2700       AliAODCaloCluster *caloCluster = new(caloClusters[jClusters++]) AliAODCaloCluster(id,
2701                                                                                         nLabel,
2702                                                                                         label,
2703                                                                                         energy,
2704                                                                                         pos,
2705                                                                                         NULL,
2706                                                                                         ttype);
2707       
2708       caloCluster->SetCaloCluster(); // to be refined!
2709
2710     } 
2711     caloClusters.Expand(jClusters); // resize TObjArray to 'remove' slots for pseudo clusters    
2712     // end of loop on calo clusters
2713
2714     // fill EMCAL cell info
2715     if (esd->GetEMCALCells()) { // protection against missing ESD information
2716       AliESDCaloCells &esdEMcells = *(esd->GetEMCALCells());
2717       Int_t nEMcell = esdEMcells.GetNumberOfCells() ;
2718       
2719       AliAODCaloCells &aodEMcells = *(aod->GetEMCALCells());
2720       aodEMcells.CreateContainer(nEMcell);
2721       aodEMcells.SetType(AliAODCaloCells::kEMCAL);
2722       for (Int_t iCell = 0; iCell < nEMcell; iCell++) {      
2723         aodEMcells.SetCell(iCell,esdEMcells.GetCellNumber(iCell),esdEMcells.GetAmplitude(iCell));
2724       }
2725       aodEMcells.Sort();
2726     }
2727
2728     // fill PHOS cell info
2729     if (esd->GetPHOSCells()) { // protection against missing ESD information
2730       AliESDCaloCells &esdPHcells = *(esd->GetPHOSCells());
2731       Int_t nPHcell = esdPHcells.GetNumberOfCells() ;
2732       
2733       AliAODCaloCells &aodPHcells = *(aod->GetPHOSCells());
2734       aodPHcells.CreateContainer(nPHcell);
2735       aodPHcells.SetType(AliAODCaloCells::kPHOS);
2736       for (Int_t iCell = 0; iCell < nPHcell; iCell++) {      
2737         aodPHcells.SetCell(iCell,esdPHcells.GetCellNumber(iCell),esdPHcells.GetAmplitude(iCell));
2738       }
2739       aodPHcells.Sort();
2740     }
2741
2742     // tracklets    
2743     AliAODTracklets &SPDTracklets = *(aod->GetTracklets());
2744     const AliMultiplicity *mult = esd->GetMultiplicity();
2745     if (mult) {
2746       if (mult->GetNumberOfTracklets()>0) {
2747         SPDTracklets.CreateContainer(mult->GetNumberOfTracklets());
2748
2749         for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
2750           SPDTracklets.SetTracklet(n, mult->GetTheta(n), mult->GetPhi(n), mult->GetDeltaPhi(n), mult->GetLabel(n));
2751         }
2752       }
2753     } else {
2754       Printf("ERROR: AliMultiplicity could not be retrieved from ESD");
2755     }
2756
2757     delete [] usedTrack;
2758     delete [] usedV0;
2759     delete [] usedKink;
2760
2761     // fill the tree for this event
2762     aodTree->Fill();
2763   } // end of event loop
2764
2765   aodTree->GetUserInfo()->Add(aod);
2766
2767   // write the tree to the specified file
2768   aodFile = aodTree->GetCurrentFile();
2769   aodFile->cd();
2770   aodTree->Write();
2771
2772   return;
2773 }
2774
2775 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2776 {
2777   // Write space-points which are then used in the alignment procedures
2778   // For the moment only ITS, TRD and TPC
2779
2780   // Load TOF clusters
2781   if (fTracker[3]){
2782     fLoader[3]->LoadRecPoints("read");
2783     TTree* tree = fLoader[3]->TreeR();
2784     if (!tree) {
2785       AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2786       return;
2787     }
2788     fTracker[3]->LoadClusters(tree);
2789   }
2790   Int_t ntracks = esd->GetNumberOfTracks();
2791   for (Int_t itrack = 0; itrack < ntracks; itrack++)
2792     {
2793       AliESDtrack *track = esd->GetTrack(itrack);
2794       Int_t nsp = 0;
2795       Int_t idx[200];
2796       for (Int_t iDet = 3; iDet >= 0; iDet--)
2797         nsp += track->GetNcls(iDet);
2798       if (nsp) {
2799         AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2800         track->SetTrackPointArray(sp);
2801         Int_t isptrack = 0;
2802         for (Int_t iDet = 3; iDet >= 0; iDet--) {
2803           AliTracker *tracker = fTracker[iDet];
2804           if (!tracker) continue;
2805           Int_t nspdet = track->GetNcls(iDet);
2806           if (nspdet <= 0) continue;
2807           track->GetClusters(iDet,idx);
2808           AliTrackPoint p;
2809           Int_t isp = 0;
2810           Int_t isp2 = 0;
2811           while (isp < nspdet) {
2812             Bool_t isvalid;
2813             TString dets = fgkDetectorName[iDet];
2814             if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2815             fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2816             fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2817             fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2818               isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2819             } else {
2820               isvalid = tracker->GetTrackPoint(idx[isp2],p); 
2821             } 
2822             isp2++;
2823             const Int_t kNTPCmax = 159;
2824             if (iDet==1 && isp2>kNTPCmax) break;   // to be fixed
2825             if (!isvalid) continue;
2826             sp->AddPoint(isptrack,&p); isptrack++; isp++;
2827           }
2828         }       
2829       }
2830     }
2831   if (fTracker[3]){
2832     fTracker[3]->UnloadClusters();
2833     fLoader[3]->UnloadRecPoints();
2834   }
2835 }
2836
2837 //_____________________________________________________________________________
2838 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2839 {
2840   // The method reads the raw-data error log
2841   // accumulated within the rawReader.
2842   // It extracts the raw-data errors related to
2843   // the current event and stores them into
2844   // a TClonesArray inside the esd object.
2845
2846   if (!fRawReader) return;
2847
2848   for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2849
2850     AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2851     if (!log) continue;
2852     if (iEvent != log->GetEventNumber()) continue;
2853
2854     esd->AddRawDataErrorLog(log);
2855   }
2856
2857 }
2858
2859 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2860   // Dump a file content into a char in TNamed
2861   ifstream in;
2862   in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2863   Int_t kBytes = (Int_t)in.tellg();
2864   printf("Size: %d \n",kBytes);
2865   TNamed *fn = 0;
2866   if(in.good()){
2867     char* memblock = new char [kBytes];
2868     in.seekg (0, ios::beg);
2869     in.read (memblock, kBytes);
2870     in.close();
2871     TString fData(memblock,kBytes);
2872     fn = new TNamed(fName,fData);
2873     printf("fData Size: %d \n",fData.Sizeof());
2874     printf("fName Size: %d \n",fName.Sizeof());
2875     printf("fn    Size: %d \n",fn->Sizeof());
2876     delete[] memblock;
2877   }
2878   else{
2879     AliInfo(Form("Could not Open %s\n",fPath.Data()));
2880   }
2881
2882   return fn;
2883 }
2884
2885 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2886   // This is not really needed in AliReconstruction at the moment
2887   // but can serve as a template
2888
2889   TList *fList = fTree->GetUserInfo();
2890   TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2891   printf("fn Size: %d \n",fn->Sizeof());
2892
2893   TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2894   const char* cdata = fn->GetTitle();
2895   printf("fTmp Size %d\n",fTmp.Sizeof());
2896
2897   int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2898   printf("calculated size %d\n",size);
2899   ofstream out(fName.Data(),ios::out | ios::binary);
2900   out.write(cdata,size);
2901   out.close();
2902
2903 }
2904   
2905 //_____________________________________________________________________________
2906 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2907 {
2908  // get the quality assurance data maker object and the loader for a detector
2909
2910   if (fQADataMaker[iDet]) 
2911     return fQADataMaker[iDet];
2912
2913   AliQADataMakerRec * qadm = NULL;
2914   if (iDet == fgkNDetectors) { //Global QA
2915      qadm = new AliGlobalQADataMaker();
2916      fQADataMaker[iDet] = qadm;
2917      return qadm;
2918   }
2919
2920   // load the QA data maker object
2921   TPluginManager* pluginManager = gROOT->GetPluginManager();
2922   TString detName = fgkDetectorName[iDet];
2923   TString qadmName = "Ali" + detName + "QADataMakerRec";
2924   if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) 
2925     return NULL;
2926
2927   // first check if a plugin is defined for the quality assurance data maker
2928   TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2929   // if not, add a plugin for it
2930   if (!pluginHandler) {
2931     AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2932     TString libs = gSystem->GetLibraries();
2933     if (libs.Contains("lib" + detName + "base.so") ||
2934         (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2935       pluginManager->AddHandler("AliQADataMakerRec", detName, 
2936                                 qadmName, detName + "qadm", qadmName + "()");
2937     } else {
2938       pluginManager->AddHandler("AliQADataMakerRec", detName, 
2939                                 qadmName, detName, qadmName + "()");
2940     }
2941     pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2942   }
2943   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2944     qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2945   }
2946
2947   fQADataMaker[iDet] = qadm;
2948
2949   return qadm;
2950 }
2951
2952 //_____________________________________________________________________________
2953 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
2954 {
2955   // run the Quality Assurance data producer
2956
2957   AliCodeTimerAuto("")
2958   TString detStr = detectors;
2959   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2960    if (!IsSelected(fgkDetectorName[iDet], detStr)) 
2961      continue;
2962    AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2963    if (!qadm) 
2964      continue;
2965    AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2966    AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2967     
2968    qadm->Exec(AliQA::kESDS, esd) ; 
2969    qadm->Increment() ; 
2970
2971    AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2972  }
2973  if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2974    AliError(Form("the following detectors were not found: %s",
2975                  detStr.Data()));
2976    if (fStopOnError) 
2977      return kFALSE;
2978  }
2979  
2980  return kTRUE;
2981   
2982 }
2983
2984 //_____________________________________________________________________________
2985 void AliReconstruction::CheckQA()
2986 {
2987 // check the QA of SIM for this run and remove the detectors 
2988 // with status Fatal
2989   
2990         TString newRunLocalReconstruction ; 
2991         TString newRunTracking ;
2992         TString newFillESD ;
2993          
2994         for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2995                 TString detName(AliQA::GetDetName(iDet)) ;
2996                 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX(iDet)) ; 
2997                 if ( qa->IsSet(AliQA::DETECTORINDEX(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2998                                 AliInfo(Form("QA status for %s in Hits and/or SDIGITS  and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2999                 } else {
3000                         if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) || 
3001                                         fRunLocalReconstruction.Contains("ALL") )  {
3002                                 newRunLocalReconstruction += detName ; 
3003                                 newRunLocalReconstruction += " " ;                      
3004                         }
3005                         if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) || 
3006                                         fRunTracking.Contains("ALL") )  {
3007                                 newRunTracking += detName ; 
3008                                 newRunTracking += " " ;                         
3009                         }
3010                         if ( fFillESD.Contains(AliQA::GetDetName(iDet)) || 
3011                                         fFillESD.Contains("ALL") )  {
3012                                 newFillESD += detName ; 
3013                                 newFillESD += " " ;                     
3014                         }
3015                 }
3016         }
3017         fRunLocalReconstruction = newRunLocalReconstruction ; 
3018         fRunTracking            = newRunTracking ; 
3019         fFillESD                = newFillESD ; 
3020 }
3021
3022 //_____________________________________________________________________________
3023 Int_t AliReconstruction::GetDetIndex(const char* detector)
3024 {
3025   // return the detector index corresponding to detector
3026   Int_t index = -1 ; 
3027   for (index = 0; index < fgkNDetectors ; index++) {
3028     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3029         break ; 
3030   }     
3031   return index ; 
3032 }