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