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