]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliReconstruction.cxx
New implementation of Print (Alberto)
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.cxx
CommitLineData
596a855f 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. //
c71de921 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// //
b26c3770 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// //
596a855f 50// The name of the galice file can be changed from the default //
e583c30d 51// "galice.root" by passing it as argument to the AliReconstruction //
52// constructor or by //
596a855f 53// //
54// rec.SetGAliceFile("..."); //
55// //
59697224 56// The local reconstruction can be switched on or off for individual //
57// detectors by //
596a855f 58// //
59697224 59// rec.SetRunLocalReconstruction("..."); //
596a855f 60// //
61// The argument is a (case sensitive) string with the names of the //
62// detectors separated by a space. The special string "ALL" selects all //
63// available detectors. This is the default. //
64// //
c71de921 65// The reconstruction of the primary vertex position can be switched off by //
66// //
67// rec.SetRunVertexFinder(kFALSE); //
68// //
b8cd5251 69// The tracking and the creation of ESD tracks can be switched on for //
70// selected detectors by //
596a855f 71// //
b8cd5251 72// rec.SetRunTracking("..."); //
596a855f 73// //
c84a5e9e 74// Uniform/nonuniform field tracking switches (default: uniform field) //
75// //
1d99986f 76// rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
c84a5e9e 77// //
596a855f 78// The filling of additional ESD information can be steered by //
79// //
80// rec.SetFillESD("..."); //
81// //
b8cd5251 82// Again, for both methods the string specifies the list of detectors. //
83// The default is "ALL". //
84// //
85// The call of the shortcut method //
86// //
87// rec.SetRunReconstruction("..."); //
88// //
89// is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
90// SetFillESD with the same detector selecting string as argument. //
596a855f 91// //
c71de921 92// The reconstruction requires digits or raw data as input. For the creation //
93// of digits and raw data have a look at the class AliSimulation. //
596a855f 94// //
24f7a148 95// For debug purposes the method SetCheckPointLevel can be used. If the //
96// argument is greater than 0, files with ESD events will be written after //
97// selected steps of the reconstruction for each event: //
98// level 1: after tracking and after filling of ESD (final) //
99// level 2: in addition after each tracking step //
100// level 3: in addition after the filling of ESD for each detector //
101// If a final check point file exists for an event, this event will be //
102// skipped in the reconstruction. The tracking and the filling of ESD for //
103// a detector will be skipped as well, if the corresponding check point //
104// file exists. The ESD event will then be loaded from the file instead. //
105// //
596a855f 106///////////////////////////////////////////////////////////////////////////////
107
024a7e64 108#include <TArrayF.h>
109#include <TFile.h>
110#include <TSystem.h>
111#include <TROOT.h>
112#include <TPluginManager.h>
fd46e2d2 113#include <TStopwatch.h>
3103d196 114#include <TGeoManager.h>
2bdb9d38 115#include <TLorentzVector.h>
596a855f 116
117#include "AliReconstruction.h"
b8cd5251 118#include "AliReconstructor.h"
815c2b38 119#include "AliLog.h"
596a855f 120#include "AliRunLoader.h"
121#include "AliRun.h"
b649205a 122#include "AliRawReaderFile.h"
123#include "AliRawReaderDate.h"
124#include "AliRawReaderRoot.h"
596a855f 125#include "AliESD.h"
1d99986f 126#include "AliESDfriend.h"
2257f27e 127#include "AliESDVertex.h"
32e449be 128#include "AliMultiplicity.h"
c84a5e9e 129#include "AliTracker.h"
2257f27e 130#include "AliVertexer.h"
c5e3e5d1 131#include "AliVertexerTracks.h"
596a855f 132#include "AliHeader.h"
133#include "AliGenEventHeader.h"
b26c3770 134#include "AliPID.h"
596a855f 135#include "AliESDpid.h"
ff8bb5ae 136#include "AliESDtrack.h"
f3a97c86 137
138#include "AliRunTag.h"
f3a97c86 139#include "AliDetectorTag.h"
140#include "AliEventTag.h"
141
98937d93 142#include "AliTrackPointArray.h"
b0314964 143#include "AliCDBManager.h"
6bae477a 144#include "AliCDBEntry.h"
145#include "AliAlignObj.h"
f3a97c86 146
b647652d 147#include "AliCentralTrigger.h"
148#include "AliCTPRawStream.h"
149
596a855f 150ClassImp(AliReconstruction)
151
152
c757bafd 153//_____________________________________________________________________________
482070f2 154const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
c757bafd 155
596a855f 156//_____________________________________________________________________________
024cf675 157AliReconstruction::AliReconstruction(const char* gAliceFilename, const char* cdbUri,
e583c30d 158 const char* name, const char* title) :
159 TNamed(name, title),
160
c84a5e9e 161 fUniformField(kTRUE),
2257f27e 162 fRunVertexFinder(kTRUE),
1f46a9ae 163 fRunHLTTracking(kFALSE),
1d99986f 164 fStopOnError(kFALSE),
165 fWriteAlignmentData(kFALSE),
166 fWriteESDfriend(kFALSE),
b647652d 167 fFillTriggerESD(kTRUE),
1d99986f 168
169 fRunLocalReconstruction("ALL"),
b8cd5251 170 fRunTracking("ALL"),
e583c30d 171 fFillESD("ALL"),
172 fGAliceFileName(gAliceFilename),
b649205a 173 fInput(""),
35042093 174 fEquipIdMap(""),
b26c3770 175 fFirstEvent(0),
176 fLastEvent(-1),
24f7a148 177 fCheckPointLevel(0),
b8cd5251 178 fOptions(),
6bae477a 179 fLoadAlignFromCDB(kTRUE),
180 fLoadAlignData("ALL"),
e583c30d 181
182 fRunLoader(NULL),
b649205a 183 fRawReader(NULL),
b8cd5251 184
98937d93 185 fVertexer(NULL),
186
6bae477a 187 fAlignObjArray(NULL),
ec92bee0 188 fCDBUri(cdbUri),
189 fSpecCDBUri()
596a855f 190{
191// create reconstruction object with default parameters
b8cd5251 192
193 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
194 fReconstructor[iDet] = NULL;
195 fLoader[iDet] = NULL;
196 fTracker[iDet] = NULL;
197 }
e47c4c2e 198 AliPID pid;
596a855f 199}
200
201//_____________________________________________________________________________
202AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
e583c30d 203 TNamed(rec),
204
c84a5e9e 205 fUniformField(rec.fUniformField),
2257f27e 206 fRunVertexFinder(rec.fRunVertexFinder),
1f46a9ae 207 fRunHLTTracking(rec.fRunHLTTracking),
1d99986f 208 fStopOnError(rec.fStopOnError),
209 fWriteAlignmentData(rec.fWriteAlignmentData),
210 fWriteESDfriend(rec.fWriteESDfriend),
b647652d 211 fFillTriggerESD(rec.fFillTriggerESD),
1d99986f 212
213 fRunLocalReconstruction(rec.fRunLocalReconstruction),
e583c30d 214 fRunTracking(rec.fRunTracking),
215 fFillESD(rec.fFillESD),
216 fGAliceFileName(rec.fGAliceFileName),
b649205a 217 fInput(rec.fInput),
35042093 218 fEquipIdMap(rec.fEquipIdMap),
b26c3770 219 fFirstEvent(rec.fFirstEvent),
220 fLastEvent(rec.fLastEvent),
24f7a148 221 fCheckPointLevel(0),
b8cd5251 222 fOptions(),
6bae477a 223 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
224 fLoadAlignData(rec.fLoadAlignData),
e583c30d 225
226 fRunLoader(NULL),
b649205a 227 fRawReader(NULL),
b8cd5251 228
98937d93 229 fVertexer(NULL),
230
6bae477a 231 fAlignObjArray(rec.fAlignObjArray),
ec92bee0 232 fCDBUri(rec.fCDBUri),
233 fSpecCDBUri()
596a855f 234{
235// copy constructor
236
ec92bee0 237 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
efd2085e 238 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
239 }
b8cd5251 240 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
241 fReconstructor[iDet] = NULL;
242 fLoader[iDet] = NULL;
243 fTracker[iDet] = NULL;
244 }
ec92bee0 245 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
246 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
247 }
596a855f 248}
249
250//_____________________________________________________________________________
251AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
252{
253// assignment operator
254
255 this->~AliReconstruction();
256 new(this) AliReconstruction(rec);
257 return *this;
258}
259
260//_____________________________________________________________________________
261AliReconstruction::~AliReconstruction()
262{
263// clean up
264
e583c30d 265 CleanUp();
efd2085e 266 fOptions.Delete();
ec92bee0 267 fSpecCDBUri.Delete();
596a855f 268}
269
024cf675 270//_____________________________________________________________________________
271void AliReconstruction::InitCDBStorage()
272{
273// activate a default CDB storage
274// First check if we have any CDB storage set, because it is used
275// to retrieve the calibration and alignment constants
276
277 AliCDBManager* man = AliCDBManager::Instance();
ec92bee0 278 if (man->IsDefaultStorageSet())
024cf675 279 {
280 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
ec92bee0 281 AliWarning("Default CDB storage has been already set !");
282 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
024cf675 283 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
ec92bee0 284 fCDBUri = "";
285 }
286 else {
287 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
288 AliWarning(Form("Default CDB storage is set to: %s",fCDBUri.Data()));
289 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
290 man->SetDefaultStorage(fCDBUri);
291 }
292
293 // Now activate the detector specific CDB storage locations
294 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
295 TString detName = fgkDetectorName[iDet];
296 TObject* obj = fSpecCDBUri.FindObject(detName.Data());
bb3da222 297 if (obj) {
298 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
299 AliWarning(Form("Specific CDB storage for %s is set to: %s",detName.Data(),obj->GetTitle()));
300 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
301 man->SetSpecificStorage(detName.Data(),obj->GetTitle());
302 }
ec92bee0 303 }
024cf675 304
305}
306
307//_____________________________________________________________________________
308void AliReconstruction::SetDefaultStorage(const char* uri) {
ec92bee0 309// Store the desired default CDB storage location
310// Activate it later within the Run() method
024cf675 311
ec92bee0 312 fCDBUri = uri;
024cf675 313
314}
315
316//_____________________________________________________________________________
317void AliReconstruction::SetSpecificStorage(const char* detName, const char* uri) {
ec92bee0 318// Store a detector-specific CDB storage location
319// Activate it later within the Run() method
024cf675 320
ec92bee0 321 TObject* obj = fSpecCDBUri.FindObject(detName);
322 if (obj) fSpecCDBUri.Remove(obj);
323 fSpecCDBUri.Add(new TNamed(detName, uri));
024cf675 324
325}
326
6bae477a 327//_____________________________________________________________________________
328Bool_t AliReconstruction::SetRunNumber()
329{
330 // The method is called in Run() in order
331 // to set a correct run number.
332 // In case of raw data reconstruction the
333 // run number is taken from the raw data header
334
335 if(AliCDBManager::Instance()->GetRun() < 0) {
336 if (!fRunLoader) {
337 AliError("No run loader is found !");
338 return kFALSE;
339 }
340 // read run number from gAlice
ec92bee0 341 if(fRunLoader->GetAliRun())
342 AliCDBManager::Instance()->SetRun(fRunLoader->GetAliRun()->GetRunNumber());
343 else {
344 if(fRawReader) {
345 if(fRawReader->NextEvent()) {
346 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
347 fRawReader->RewindEvents();
348 }
349 else {
350 AliError("No raw-data events found !");
351 return kFALSE;
352 }
353 }
354 else {
355 AliError("Neither gAlice nor RawReader objects are found !");
356 return kFALSE;
357 }
358 }
359 AliInfo(Form("CDB Run number: %d",AliCDBManager::Instance()->GetRun()));
6bae477a 360 }
361 return kTRUE;
362}
363
364//_____________________________________________________________________________
365Bool_t AliReconstruction::ApplyAlignObjsToGeom(TObjArray* alObjArray)
366{
367 // Read collection of alignment objects (AliAlignObj derived) saved
368 // in the TClonesArray ClArrayName and apply them to the geometry
369 // manager singleton.
370 //
371 alObjArray->Sort();
372 Int_t nvols = alObjArray->GetEntriesFast();
373
dc0984f8 374 Bool_t flag = kTRUE;
375
6bae477a 376 for(Int_t j=0; j<nvols; j++)
377 {
378 AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
dc0984f8 379 if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
6bae477a 380 }
381
382 if (AliDebugLevelClass() >= 1) {
e0625db4 383 gGeoManager->GetTopNode()->CheckOverlaps(20);
6bae477a 384 TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
385 if(ovexlist->GetEntriesFast()){
386 AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
387 }
388 }
389
dc0984f8 390 return flag;
6bae477a 391
392}
393
394//_____________________________________________________________________________
395Bool_t AliReconstruction::SetAlignObjArraySingleDet(const char* detName)
396{
397 // Fills array of single detector's alignable objects from CDB
398
399 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
400
401 AliCDBEntry *entry;
402
403 AliCDBPath path(detName,"Align","Data");
404
405 entry=AliCDBManager::Instance()->Get(path.GetPath());
406 if(!entry){
407 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
408 return kFALSE;
409 }
410 entry->SetOwner(1);
411 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
412 alignArray->SetOwner(0);
413 AliDebug(2,Form("Found %d alignment objects for %s",
414 alignArray->GetEntries(),detName));
415
416 AliAlignObj *alignObj=0;
417 TIter iter(alignArray);
418
419 // loop over align objects in detector
420 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
421 fAlignObjArray->Add(alignObj);
422 }
423 // delete entry --- Don't delete, it is cached!
424
425 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
426 return kTRUE;
427
428}
429
430//_____________________________________________________________________________
431Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
432{
433 // Read the alignment objects from CDB.
434 // Each detector is supposed to have the
435 // alignment objects in DET/Align/Data CDB path.
436 // All the detector objects are then collected,
437 // sorted by geometry level (starting from ALIC) and
438 // then applied to the TGeo geometry.
439 // Finally an overlaps check is performed.
440
441 // Load alignment data from CDB and fill fAlignObjArray
442 if(fLoadAlignFromCDB){
443 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
444
445 //fAlignObjArray->RemoveAll();
446 fAlignObjArray->Clear();
447 fAlignObjArray->SetOwner(0);
448
449 TString detStr = detectors;
450 TString dataNotLoaded="";
451 TString dataLoaded="";
452
453 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
454 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
455 if(!SetAlignObjArraySingleDet(fgkDetectorName[iDet])){
456 dataNotLoaded += fgkDetectorName[iDet];
457 dataNotLoaded += " ";
458 } else {
459 dataLoaded += fgkDetectorName[iDet];
460 dataLoaded += " ";
461 }
462 } // end loop over detectors
463
464 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
465 dataNotLoaded += detStr;
466 AliInfo(Form("Alignment data loaded for: %s",
467 dataLoaded.Data()));
468 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
469 dataNotLoaded.Data()));
470 } // fLoadAlignFromCDB flag
471
472 // Check if the array with alignment objects was
473 // provided by the user. If yes, apply the objects
474 // to the present TGeo geometry
475 if (fAlignObjArray) {
476 if (gGeoManager && gGeoManager->IsClosed()) {
477 if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
dc0984f8 478 AliError("The misalignment of one or more volumes failed!"
479 "Compare the list of simulated detectors and the list of detector alignment data!");
6bae477a 480 return kFALSE;
481 }
482 }
483 else {
484 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
485 return kFALSE;
486 }
487 }
488
489 return kTRUE;
490}
596a855f 491
492//_____________________________________________________________________________
493void AliReconstruction::SetGAliceFile(const char* fileName)
494{
495// set the name of the galice file
496
497 fGAliceFileName = fileName;
498}
499
efd2085e 500//_____________________________________________________________________________
501void AliReconstruction::SetOption(const char* detector, const char* option)
502{
503// set options for the reconstruction of a detector
504
505 TObject* obj = fOptions.FindObject(detector);
506 if (obj) fOptions.Remove(obj);
507 fOptions.Add(new TNamed(detector, option));
508}
509
596a855f 510
511//_____________________________________________________________________________
b26c3770 512Bool_t AliReconstruction::Run(const char* input,
513 Int_t firstEvent, Int_t lastEvent)
596a855f 514{
515// run the reconstruction
516
b649205a 517 // set the input
518 if (!input) input = fInput.Data();
519 TString fileName(input);
520 if (fileName.EndsWith("/")) {
521 fRawReader = new AliRawReaderFile(fileName);
522 } else if (fileName.EndsWith(".root")) {
523 fRawReader = new AliRawReaderRoot(fileName);
524 } else if (!fileName.IsNull()) {
525 fRawReader = new AliRawReaderDate(fileName);
526 fRawReader->SelectEvents(7);
527 }
35042093 528 if (!fEquipIdMap.IsNull() && fRawReader)
529 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
530
b649205a 531
f08fc9f5 532 // get the run loader
533 if (!InitRunLoader()) return kFALSE;
596a855f 534
ec92bee0 535 // Initialize the CDB storage
536 InitCDBStorage();
537
6bae477a 538 // Set run number in CDBManager (if it is not already set by the user)
539 if (!SetRunNumber()) if (fStopOnError) return kFALSE;
540
541 // Import ideal TGeo geometry and apply misalignment
542 if (!gGeoManager) {
543 TString geom(gSystem->DirName(fGAliceFileName));
544 geom += "/geometry.root";
545 TGeoManager::Import(geom.Data());
546 if (!gGeoManager) if (fStopOnError) return kFALSE;
547 }
548 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
549
550 // Temporary fix by A.Gheata
551 // Could be removed with the next Root version (>5.11)
552 if (gGeoManager) {
553 TIter next(gGeoManager->GetListOfVolumes());
554 TGeoVolume *vol;
555 while ((vol = (TGeoVolume *)next())) {
556 if (vol->GetVoxels()) {
557 if (vol->GetVoxels()->NeedRebuild()) {
558 vol->GetVoxels()->Voxelize();
559 vol->FindOverlaps();
560 }
561 }
562 }
563 }
564
596a855f 565 // local reconstruction
59697224 566 if (!fRunLocalReconstruction.IsNull()) {
567 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
e583c30d 568 if (fStopOnError) {CleanUp(); return kFALSE;}
596a855f 569 }
570 }
b26c3770 571// if (!fRunVertexFinder && fRunTracking.IsNull() &&
572// fFillESD.IsNull()) return kTRUE;
2257f27e 573
574 // get vertexer
575 if (fRunVertexFinder && !CreateVertexer()) {
576 if (fStopOnError) {
577 CleanUp();
578 return kFALSE;
579 }
580 }
596a855f 581
f08fc9f5 582 // get trackers
b8cd5251 583 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
24f7a148 584 if (fStopOnError) {
585 CleanUp();
586 return kFALSE;
587 }
596a855f 588 }
24f7a148 589
9db3a215 590
591 TStopwatch stopwatch;
592 stopwatch.Start();
593
b26c3770 594 // get the possibly already existing ESD file and tree
1f46a9ae 595 AliESD* esd = new AliESD; AliESD* hltesd = new AliESD;
b26c3770 596 TFile* fileOld = NULL;
1f46a9ae 597 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
b26c3770 598 if (!gSystem->AccessPathName("AliESDs.root")){
599 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
600 fileOld = TFile::Open("AliESDs.old.root");
601 if (fileOld && fileOld->IsOpen()) {
602 treeOld = (TTree*) fileOld->Get("esdTree");
603 if (treeOld) treeOld->SetBranchAddress("ESD", &esd);
1f46a9ae 604 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
605 if (hlttreeOld) hlttreeOld->SetBranchAddress("ESD", &hltesd);
b26c3770 606 }
607 }
608
36711aa4 609 // create the ESD output file and tree
596a855f 610 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
611 if (!file->IsOpen()) {
815c2b38 612 AliError("opening AliESDs.root failed");
b26c3770 613 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
596a855f 614 }
36711aa4 615 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
616 tree->Branch("ESD", "AliESD", &esd);
1f46a9ae 617 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
618 hlttree->Branch("ESD", "AliESD", &hltesd);
619 delete esd; delete hltesd;
620 esd = NULL; hltesd = NULL;
1d99986f 621
622 // create the file and tree with ESD additions
623 TFile *filef=0; TTree *treef=0; AliESDfriend *esdf=0;
624 if (fWriteESDfriend) {
625 filef = TFile::Open("AliESDfriends.root", "RECREATE");
626 if (!filef->IsOpen()) {
627 AliError("opening AliESDfriends.root failed");
628 }
629 treef = new TTree("esdFriendTree", "Tree with ESD friends");
630 treef->Branch("ESDfriend", "AliESDfriend", &esdf);
631 }
632
36711aa4 633 gROOT->cd();
596a855f 634
c5e3e5d1 635 AliVertexerTracks tVertexer;
636
596a855f 637 // loop over events
b649205a 638 if (fRawReader) fRawReader->RewindEvents();
f08fc9f5 639
596a855f 640 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
b26c3770 641 if (fRawReader) fRawReader->NextEvent();
642 if ((iEvent < firstEvent) || ((lastEvent >= 0) && (iEvent > lastEvent))) {
643 // copy old ESD to the new one
644 if (treeOld) {
645 treeOld->SetBranchAddress("ESD", &esd);
646 treeOld->GetEntry(iEvent);
647 }
648 tree->Fill();
1f46a9ae 649 if (hlttreeOld) {
650 hlttreeOld->SetBranchAddress("ESD", &hltesd);
651 hlttreeOld->GetEntry(iEvent);
652 }
653 hlttree->Fill();
b26c3770 654 continue;
655 }
656
815c2b38 657 AliInfo(Form("processing event %d", iEvent));
596a855f 658 fRunLoader->GetEvent(iEvent);
24f7a148 659
660 char fileName[256];
661 sprintf(fileName, "ESD_%d.%d_final.root",
f08fc9f5 662 fRunLoader->GetHeader()->GetRun(),
663 fRunLoader->GetHeader()->GetEventNrInRun());
24f7a148 664 if (!gSystem->AccessPathName(fileName)) continue;
665
b26c3770 666 // local reconstruction
667 if (!fRunLocalReconstruction.IsNull()) {
668 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
669 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
670 }
671 }
672
1f46a9ae 673 esd = new AliESD; hltesd = new AliESD;
f08fc9f5 674 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1f46a9ae 675 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
f08fc9f5 676 esd->SetEventNumber(fRunLoader->GetHeader()->GetEventNrInRun());
1f46a9ae 677 hltesd->SetEventNumber(fRunLoader->GetHeader()->GetEventNrInRun());
d6ee376f 678
679 // Set magnetic field from the tracker
680 esd->SetMagneticField(AliTracker::GetBz());
681 hltesd->SetMagneticField(AliTracker::GetBz());
596a855f 682
2257f27e 683 // vertex finder
684 if (fRunVertexFinder) {
685 if (!ReadESD(esd, "vertex")) {
686 if (!RunVertexFinder(esd)) {
b26c3770 687 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
2257f27e 688 }
689 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
690 }
691 }
692
1f46a9ae 693 // HLT tracking
694 if (!fRunTracking.IsNull()) {
695 if (fRunHLTTracking) {
696 hltesd->SetVertex(esd->GetVertex());
697 if (!RunHLTTracking(hltesd)) {
698 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
699 }
700 }
701 }
702
596a855f 703 // barrel tracking
b8cd5251 704 if (!fRunTracking.IsNull()) {
24f7a148 705 if (!ReadESD(esd, "tracking")) {
706 if (!RunTracking(esd)) {
b26c3770 707 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
24f7a148 708 }
709 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
596a855f 710 }
711 }
712
713 // fill ESD
714 if (!fFillESD.IsNull()) {
715 if (!FillESD(esd, fFillESD)) {
b26c3770 716 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
596a855f 717 }
718 }
719
720 // combined PID
721 AliESDpid::MakePID(esd);
24f7a148 722 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
596a855f 723
b647652d 724 if (fFillTriggerESD) {
725 if (!ReadESD(esd, "trigger")) {
726 if (!FillTriggerESD(esd)) {
727 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
728 }
729 if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
730 }
731 }
732
a28195f7 733 esd->SetPrimaryVertex(tVertexer.FindPrimaryVertex(esd));
c5e3e5d1 734
596a855f 735 // write ESD
36711aa4 736 tree->Fill();
1f46a9ae 737 // write HLT ESD
738 hlttree->Fill();
24f7a148 739
1d99986f 740 // write ESD friend
741 if (fWriteESDfriend) {
d75007f6 742 esdf=new AliESDfriend();
743 esd->GetESDfriend(esdf);
1d99986f 744 treef->Fill();
745 }
746
f3a97c86 747 if (fCheckPointLevel > 0) WriteESD(esd, "final");
748
1f46a9ae 749 delete esd; delete hltesd;
750 esd = NULL; hltesd = NULL;
596a855f 751 }
752
9db3a215 753 AliInfo(Form("Execution time for filling ESD : R:%.2fs C:%.2fs",
754 stopwatch.RealTime(),stopwatch.CpuTime()));
755
36711aa4 756 file->cd();
757 tree->Write();
1f46a9ae 758 hlttree->Write();
f3a97c86 759
1d99986f 760 if (fWriteESDfriend) {
761 filef->cd();
762 treef->Write(); delete treef; filef->Close(); delete filef;
763 }
764
f3a97c86 765 // Create tags for the events in the ESD tree (the ESD tree is always present)
766 // In case of empty events the tags will contain dummy values
767 CreateTag(file);
b26c3770 768 CleanUp(file, fileOld);
596a855f 769
770 return kTRUE;
771}
772
773
774//_____________________________________________________________________________
59697224 775Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
596a855f 776{
59697224 777// run the local reconstruction
596a855f 778
030b532d 779 TStopwatch stopwatch;
780 stopwatch.Start();
781
596a855f 782 TString detStr = detectors;
b8cd5251 783 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
784 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
785 AliReconstructor* reconstructor = GetReconstructor(iDet);
786 if (!reconstructor) continue;
b26c3770 787 if (reconstructor->HasLocalReconstruction()) continue;
b8cd5251 788
789 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
790 TStopwatch stopwatchDet;
791 stopwatchDet.Start();
792 if (fRawReader) {
793 fRawReader->RewindEvents();
794 reconstructor->Reconstruct(fRunLoader, fRawReader);
795 } else {
796 reconstructor->Reconstruct(fRunLoader);
596a855f 797 }
5f8272e1 798 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
799 fgkDetectorName[iDet],
800 stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
596a855f 801 }
802
803 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
815c2b38 804 AliError(Form("the following detectors were not found: %s",
805 detStr.Data()));
596a855f 806 if (fStopOnError) return kFALSE;
807 }
808
5f8272e1 809 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
810 stopwatch.RealTime(),stopwatch.CpuTime()));
030b532d 811
596a855f 812 return kTRUE;
813}
814
b26c3770 815//_____________________________________________________________________________
816Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
817{
818// run the local reconstruction
819
820 TStopwatch stopwatch;
821 stopwatch.Start();
822
823 TString detStr = detectors;
824 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
825 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
826 AliReconstructor* reconstructor = GetReconstructor(iDet);
827 if (!reconstructor) continue;
828 AliLoader* loader = fLoader[iDet];
829
830 // conversion of digits
831 if (fRawReader && reconstructor->HasDigitConversion()) {
832 AliInfo(Form("converting raw data digits into root objects for %s",
833 fgkDetectorName[iDet]));
834 TStopwatch stopwatchDet;
835 stopwatchDet.Start();
836 loader->LoadDigits("update");
837 loader->CleanDigits();
838 loader->MakeDigitsContainer();
839 TTree* digitsTree = loader->TreeD();
840 reconstructor->ConvertDigits(fRawReader, digitsTree);
841 loader->WriteDigits("OVERWRITE");
842 loader->UnloadDigits();
5f8272e1 843 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
844 fgkDetectorName[iDet],
845 stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
b26c3770 846 }
847
848 // local reconstruction
849 if (!reconstructor->HasLocalReconstruction()) continue;
850 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
851 TStopwatch stopwatchDet;
852 stopwatchDet.Start();
853 loader->LoadRecPoints("update");
854 loader->CleanRecPoints();
855 loader->MakeRecPointsContainer();
856 TTree* clustersTree = loader->TreeR();
857 if (fRawReader && !reconstructor->HasDigitConversion()) {
858 reconstructor->Reconstruct(fRawReader, clustersTree);
859 } else {
860 loader->LoadDigits("read");
861 TTree* digitsTree = loader->TreeD();
862 if (!digitsTree) {
863 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
864 if (fStopOnError) return kFALSE;
865 } else {
866 reconstructor->Reconstruct(digitsTree, clustersTree);
867 }
868 loader->UnloadDigits();
869 }
870 loader->WriteRecPoints("OVERWRITE");
871 loader->UnloadRecPoints();
5f8272e1 872 AliDebug(1,Form("Execution time for %s: R:%.2fs C:%.2fs",
873 fgkDetectorName[iDet],
874 stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
b26c3770 875 }
876
877 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
878 AliError(Form("the following detectors were not found: %s",
879 detStr.Data()));
880 if (fStopOnError) return kFALSE;
881 }
5f8272e1 882
883 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
884 stopwatch.RealTime(),stopwatch.CpuTime()));
b26c3770 885
886 return kTRUE;
887}
888
596a855f 889//_____________________________________________________________________________
2257f27e 890Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
596a855f 891{
892// run the barrel tracking
893
030b532d 894 TStopwatch stopwatch;
895 stopwatch.Start();
896
2257f27e 897 AliESDVertex* vertex = NULL;
898 Double_t vtxPos[3] = {0, 0, 0};
899 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
900 TArrayF mcVertex(3);
a6b0b91b 901 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
902 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
903 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
904 }
2257f27e 905
b8cd5251 906 if (fVertexer) {
815c2b38 907 AliInfo("running the ITS vertex finder");
b26c3770 908 if (fLoader[0]) fLoader[0]->LoadRecPoints();
b8cd5251 909 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
b26c3770 910 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
2257f27e 911 if(!vertex){
815c2b38 912 AliWarning("Vertex not found");
c710f220 913 vertex = new AliESDVertex();
d1a50cb5 914 vertex->SetName("default");
2257f27e 915 }
916 else {
917 vertex->SetTruePos(vtxPos); // store also the vertex from MC
d1a50cb5 918 vertex->SetName("reconstructed");
2257f27e 919 }
920
921 } else {
815c2b38 922 AliInfo("getting the primary vertex from MC");
2257f27e 923 vertex = new AliESDVertex(vtxPos, vtxErr);
924 }
925
926 if (vertex) {
927 vertex->GetXYZ(vtxPos);
928 vertex->GetSigmaXYZ(vtxErr);
929 } else {
815c2b38 930 AliWarning("no vertex reconstructed");
2257f27e 931 vertex = new AliESDVertex(vtxPos, vtxErr);
932 }
933 esd->SetVertex(vertex);
32e449be 934 // if SPD multiplicity has been determined, it is stored in the ESD
935 AliMultiplicity *mult= fVertexer->GetMultiplicity();
936 if(mult)esd->SetMultiplicity(mult);
937
b8cd5251 938 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
939 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
940 }
2257f27e 941 delete vertex;
942
5f8272e1 943 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
944 stopwatch.RealTime(),stopwatch.CpuTime()));
2257f27e 945
946 return kTRUE;
947}
948
1f46a9ae 949//_____________________________________________________________________________
950Bool_t AliReconstruction::RunHLTTracking(AliESD*& esd)
951{
952// run the HLT barrel tracking
953
954 TStopwatch stopwatch;
955 stopwatch.Start();
956
957 if (!fRunLoader) {
958 AliError("Missing runLoader!");
959 return kFALSE;
960 }
961
962 AliInfo("running HLT tracking");
963
964 // Get a pointer to the HLT reconstructor
965 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
966 if (!reconstructor) return kFALSE;
967
968 // TPC + ITS
969 for (Int_t iDet = 1; iDet >= 0; iDet--) {
970 TString detName = fgkDetectorName[iDet];
971 AliDebug(1, Form("%s HLT tracking", detName.Data()));
972 reconstructor->SetOption(detName.Data());
973 AliTracker *tracker = reconstructor->CreateTracker(fRunLoader);
974 if (!tracker) {
975 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
976 if (fStopOnError) return kFALSE;
9dcc06e1 977 continue;
1f46a9ae 978 }
979 Double_t vtxPos[3];
980 Double_t vtxErr[3]={0.005,0.005,0.010};
981 const AliESDVertex *vertex = esd->GetVertex();
982 vertex->GetXYZ(vtxPos);
983 tracker->SetVertex(vtxPos,vtxErr);
984 if(iDet != 1) {
985 fLoader[iDet]->LoadRecPoints("read");
986 TTree* tree = fLoader[iDet]->TreeR();
987 if (!tree) {
988 AliError(Form("Can't get the %s cluster tree", detName.Data()));
989 return kFALSE;
990 }
991 tracker->LoadClusters(tree);
992 }
993 if (tracker->Clusters2Tracks(esd) != 0) {
994 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
995 return kFALSE;
996 }
997 if(iDet != 1) {
998 tracker->UnloadClusters();
999 }
1000 delete tracker;
1001 }
1002
5f8272e1 1003 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1004 stopwatch.RealTime(),stopwatch.CpuTime()));
1f46a9ae 1005
1006 return kTRUE;
1007}
1008
2257f27e 1009//_____________________________________________________________________________
1010Bool_t AliReconstruction::RunTracking(AliESD*& esd)
1011{
1012// run the barrel tracking
1013
1014 TStopwatch stopwatch;
1015 stopwatch.Start();
24f7a148 1016
815c2b38 1017 AliInfo("running tracking");
596a855f 1018
b8cd5251 1019 // pass 1: TPC + ITS inwards
1020 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1021 if (!fTracker[iDet]) continue;
1022 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
24f7a148 1023
b8cd5251 1024 // load clusters
1025 fLoader[iDet]->LoadRecPoints("read");
1026 TTree* tree = fLoader[iDet]->TreeR();
1027 if (!tree) {
1028 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
24f7a148 1029 return kFALSE;
1030 }
b8cd5251 1031 fTracker[iDet]->LoadClusters(tree);
1032
1033 // run tracking
1034 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1035 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
24f7a148 1036 return kFALSE;
1037 }
b8cd5251 1038 if (fCheckPointLevel > 1) {
1039 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1040 }
878e1fe1 1041 // preliminary PID in TPC needed by the ITS tracker
1042 if (iDet == 1) {
1043 GetReconstructor(1)->FillESD(fRunLoader, esd);
b26c3770 1044 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
878e1fe1 1045 AliESDpid::MakePID(esd);
1046 }
b8cd5251 1047 }
596a855f 1048
b8cd5251 1049 // pass 2: ALL backwards
1050 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1051 if (!fTracker[iDet]) continue;
1052 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1053
1054 // load clusters
1055 if (iDet > 1) { // all except ITS, TPC
1056 TTree* tree = NULL;
7b61cd9c 1057 fLoader[iDet]->LoadRecPoints("read");
1058 tree = fLoader[iDet]->TreeR();
b8cd5251 1059 if (!tree) {
1060 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
24f7a148 1061 return kFALSE;
1062 }
b8cd5251 1063 fTracker[iDet]->LoadClusters(tree);
1064 }
24f7a148 1065
b8cd5251 1066 // run tracking
1067 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1068 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1069 return kFALSE;
1070 }
1071 if (fCheckPointLevel > 1) {
1072 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1073 }
24f7a148 1074
b8cd5251 1075 // unload clusters
1076 if (iDet > 2) { // all except ITS, TPC, TRD
1077 fTracker[iDet]->UnloadClusters();
7b61cd9c 1078 fLoader[iDet]->UnloadRecPoints();
b8cd5251 1079 }
8f37df88 1080 // updated PID in TPC needed by the ITS tracker -MI
1081 if (iDet == 1) {
1082 GetReconstructor(1)->FillESD(fRunLoader, esd);
1083 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1084 AliESDpid::MakePID(esd);
1085 }
b8cd5251 1086 }
596a855f 1087
98937d93 1088 // write space-points to the ESD in case alignment data output
1089 // is switched on
1090 if (fWriteAlignmentData)
1091 WriteAlignmentData(esd);
1092
b8cd5251 1093 // pass 3: TRD + TPC + ITS refit inwards
1094 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1095 if (!fTracker[iDet]) continue;
1096 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
596a855f 1097
b8cd5251 1098 // run tracking
1099 if (fTracker[iDet]->RefitInward(esd) != 0) {
1100 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1101 return kFALSE;
1102 }
1103 if (fCheckPointLevel > 1) {
1104 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1105 }
596a855f 1106
b8cd5251 1107 // unload clusters
1108 fTracker[iDet]->UnloadClusters();
1109 fLoader[iDet]->UnloadRecPoints();
1110 }
ff8bb5ae 1111 //
1112 // Propagate track to the vertex - if not done by ITS
1113 //
1114 Int_t ntracks = esd->GetNumberOfTracks();
1115 for (Int_t itrack=0; itrack<ntracks; itrack++){
1116 const Double_t kRadius = 3; // beam pipe radius
1117 const Double_t kMaxStep = 5; // max step
1118 const Double_t kMaxD = 123456; // max distance to prim vertex
1119 Double_t fieldZ = AliTracker::GetBz(); //
1120 AliESDtrack * track = esd->GetTrack(itrack);
1121 if (!track) continue;
1122 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
a28195f7 1123 track->PropagateTo(kRadius, fieldZ, track->GetMass(),kMaxStep,kTRUE);
ff8bb5ae 1124 track->RelateToVertex(esd->GetVertex(),fieldZ, kMaxD);
1125 }
1126
5f8272e1 1127 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1128 stopwatch.RealTime(),stopwatch.CpuTime()));
030b532d 1129
596a855f 1130 return kTRUE;
1131}
1132
1133//_____________________________________________________________________________
24f7a148 1134Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
596a855f 1135{
1136// fill the event summary data
1137
030b532d 1138 TStopwatch stopwatch;
1139 stopwatch.Start();
815c2b38 1140 AliInfo("filling ESD");
030b532d 1141
596a855f 1142 TString detStr = detectors;
b8cd5251 1143 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1144 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1145 AliReconstructor* reconstructor = GetReconstructor(iDet);
1146 if (!reconstructor) continue;
1147
1148 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1149 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
b26c3770 1150 TTree* clustersTree = NULL;
1151 if (reconstructor->HasLocalReconstruction() && fLoader[iDet]) {
1152 fLoader[iDet]->LoadRecPoints("read");
1153 clustersTree = fLoader[iDet]->TreeR();
1154 if (!clustersTree) {
1155 AliError(Form("Can't get the %s clusters tree",
1156 fgkDetectorName[iDet]));
1157 if (fStopOnError) return kFALSE;
1158 }
1159 }
1160 if (fRawReader && !reconstructor->HasDigitConversion()) {
1161 reconstructor->FillESD(fRawReader, clustersTree, esd);
1162 } else {
1163 TTree* digitsTree = NULL;
1164 if (fLoader[iDet]) {
1165 fLoader[iDet]->LoadDigits("read");
1166 digitsTree = fLoader[iDet]->TreeD();
1167 if (!digitsTree) {
1168 AliError(Form("Can't get the %s digits tree",
1169 fgkDetectorName[iDet]));
1170 if (fStopOnError) return kFALSE;
1171 }
1172 }
1173 reconstructor->FillESD(digitsTree, clustersTree, esd);
1174 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1175 }
1176 if (reconstructor->HasLocalReconstruction() && fLoader[iDet]) {
1177 fLoader[iDet]->UnloadRecPoints();
1178 }
1179
b8cd5251 1180 if (fRawReader) {
1181 reconstructor->FillESD(fRunLoader, fRawReader, esd);
1182 } else {
1183 reconstructor->FillESD(fRunLoader, esd);
24f7a148 1184 }
b8cd5251 1185 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
596a855f 1186 }
1187 }
1188
1189 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
815c2b38 1190 AliError(Form("the following detectors were not found: %s",
1191 detStr.Data()));
596a855f 1192 if (fStopOnError) return kFALSE;
1193 }
1194
5f8272e1 1195 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1196 stopwatch.RealTime(),stopwatch.CpuTime()));
030b532d 1197
596a855f 1198 return kTRUE;
1199}
1200
b647652d 1201//_____________________________________________________________________________
1202Bool_t AliReconstruction::FillTriggerESD(AliESD*& esd)
1203{
1204 // Reads the trigger decision which is
1205 // stored in Trigger.root file and fills
1206 // the corresponding esd entries
1207
1208 AliInfo("Filling trigger information into the ESD");
1209
1210 if (fRawReader) {
1211 AliCTPRawStream input(fRawReader);
1212 if (!input.Next()) {
1213 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1214 return kFALSE;
1215 }
1216 esd->SetTriggerMask(input.GetClassMask());
1217 esd->SetTriggerCluster(input.GetClusterMask());
1218 }
1219 else {
1220 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1221 if (runloader) {
1222 if (!runloader->LoadTrigger()) {
1223 AliCentralTrigger *aCTP = runloader->GetTrigger();
1224 esd->SetTriggerMask(aCTP->GetClassMask());
1225 esd->SetTriggerCluster(aCTP->GetClusterMask());
1226 }
1227 else {
1228 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1229 return kFALSE;
1230 }
1231 }
1232 else {
1233 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1234 return kFALSE;
1235 }
1236 }
1237
1238 return kTRUE;
1239}
596a855f 1240
1241//_____________________________________________________________________________
1242Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1243{
1244// check whether detName is contained in detectors
1245// if yes, it is removed from detectors
1246
1247 // check if all detectors are selected
1248 if ((detectors.CompareTo("ALL") == 0) ||
1249 detectors.BeginsWith("ALL ") ||
1250 detectors.EndsWith(" ALL") ||
1251 detectors.Contains(" ALL ")) {
1252 detectors = "ALL";
1253 return kTRUE;
1254 }
1255
1256 // search for the given detector
1257 Bool_t result = kFALSE;
1258 if ((detectors.CompareTo(detName) == 0) ||
1259 detectors.BeginsWith(detName+" ") ||
1260 detectors.EndsWith(" "+detName) ||
1261 detectors.Contains(" "+detName+" ")) {
1262 detectors.ReplaceAll(detName, "");
1263 result = kTRUE;
1264 }
1265
1266 // clean up the detectors string
1267 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1268 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1269 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1270
1271 return result;
1272}
e583c30d 1273
f08fc9f5 1274//_____________________________________________________________________________
1275Bool_t AliReconstruction::InitRunLoader()
1276{
1277// get or create the run loader
1278
1279 if (gAlice) delete gAlice;
1280 gAlice = NULL;
1281
b26c3770 1282 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1283 // load all base libraries to get the loader classes
1284 TString libs = gSystem->GetLibraries();
1285 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1286 TString detName = fgkDetectorName[iDet];
1287 if (detName == "HLT") continue;
1288 if (libs.Contains("lib" + detName + "base.so")) continue;
1289 gSystem->Load("lib" + detName + "base.so");
1290 }
f08fc9f5 1291 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1292 if (!fRunLoader) {
1293 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1294 CleanUp();
1295 return kFALSE;
1296 }
b26c3770 1297 fRunLoader->CdGAFile();
1298 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1299 if (fRunLoader->LoadgAlice() == 0) {
1300 gAlice = fRunLoader->GetAliRun();
c84a5e9e 1301 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
b26c3770 1302 }
f08fc9f5 1303 }
1304 if (!gAlice && !fRawReader) {
1305 AliError(Form("no gAlice object found in file %s",
1306 fGAliceFileName.Data()));
1307 CleanUp();
1308 return kFALSE;
1309 }
1310
1311 } else { // galice.root does not exist
1312 if (!fRawReader) {
1313 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1314 CleanUp();
1315 return kFALSE;
1316 }
1317 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1318 AliConfig::GetDefaultEventFolderName(),
1319 "recreate");
1320 if (!fRunLoader) {
1321 AliError(Form("could not create run loader in file %s",
1322 fGAliceFileName.Data()));
1323 CleanUp();
1324 return kFALSE;
1325 }
1326 fRunLoader->MakeTree("E");
1327 Int_t iEvent = 0;
1328 while (fRawReader->NextEvent()) {
1329 fRunLoader->SetEventNumber(iEvent);
1330 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1331 iEvent, iEvent);
1332 fRunLoader->MakeTree("H");
1333 fRunLoader->TreeE()->Fill();
1334 iEvent++;
1335 }
1336 fRawReader->RewindEvents();
1337 fRunLoader->WriteHeader("OVERWRITE");
1338 fRunLoader->CdGAFile();
1339 fRunLoader->Write(0, TObject::kOverwrite);
1340// AliTracker::SetFieldMap(???);
1341 }
1342
1343 return kTRUE;
1344}
1345
c757bafd 1346//_____________________________________________________________________________
b8cd5251 1347AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
c757bafd 1348{
f08fc9f5 1349// get the reconstructor object and the loader for a detector
c757bafd 1350
b8cd5251 1351 if (fReconstructor[iDet]) return fReconstructor[iDet];
1352
1353 // load the reconstructor object
1354 TPluginManager* pluginManager = gROOT->GetPluginManager();
1355 TString detName = fgkDetectorName[iDet];
1356 TString recName = "Ali" + detName + "Reconstructor";
f08fc9f5 1357 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
b8cd5251 1358
1359 if (detName == "HLT") {
1360 if (!gROOT->GetClass("AliLevel3")) {
1361 gSystem->Load("libAliL3Src.so");
1362 gSystem->Load("libAliL3Misc.so");
1363 gSystem->Load("libAliL3Hough.so");
1364 gSystem->Load("libAliL3Comp.so");
1365 }
1366 }
1367
1368 AliReconstructor* reconstructor = NULL;
1369 // first check if a plugin is defined for the reconstructor
1370 TPluginHandler* pluginHandler =
1371 pluginManager->FindHandler("AliReconstructor", detName);
f08fc9f5 1372 // if not, add a plugin for it
1373 if (!pluginHandler) {
b8cd5251 1374 AliDebug(1, Form("defining plugin for %s", recName.Data()));
b26c3770 1375 TString libs = gSystem->GetLibraries();
1376 if (libs.Contains("lib" + detName + "base.so") ||
1377 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
b8cd5251 1378 pluginManager->AddHandler("AliReconstructor", detName,
1379 recName, detName + "rec", recName + "()");
1380 } else {
1381 pluginManager->AddHandler("AliReconstructor", detName,
1382 recName, detName, recName + "()");
c757bafd 1383 }
b8cd5251 1384 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1385 }
1386 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1387 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
c757bafd 1388 }
b8cd5251 1389 if (reconstructor) {
1390 TObject* obj = fOptions.FindObject(detName.Data());
1391 if (obj) reconstructor->SetOption(obj->GetTitle());
b26c3770 1392 reconstructor->Init(fRunLoader);
b8cd5251 1393 fReconstructor[iDet] = reconstructor;
1394 }
1395
f08fc9f5 1396 // get or create the loader
1397 if (detName != "HLT") {
1398 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1399 if (!fLoader[iDet]) {
1400 AliConfig::Instance()
1401 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1402 detName, detName);
1403 // first check if a plugin is defined for the loader
1404 TPluginHandler* pluginHandler =
1405 pluginManager->FindHandler("AliLoader", detName);
1406 // if not, add a plugin for it
1407 if (!pluginHandler) {
1408 TString loaderName = "Ali" + detName + "Loader";
1409 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1410 pluginManager->AddHandler("AliLoader", detName,
1411 loaderName, detName + "base",
1412 loaderName + "(const char*, TFolder*)");
1413 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1414 }
1415 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1416 fLoader[iDet] =
1417 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1418 fRunLoader->GetEventFolder());
1419 }
1420 if (!fLoader[iDet]) { // use default loader
1421 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1422 }
1423 if (!fLoader[iDet]) {
1424 AliWarning(Form("couldn't get loader for %s", detName.Data()));
6667b602 1425 if (fStopOnError) return NULL;
f08fc9f5 1426 } else {
1427 fRunLoader->AddLoader(fLoader[iDet]);
1428 fRunLoader->CdGAFile();
1429 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1430 fRunLoader->Write(0, TObject::kOverwrite);
1431 }
1432 }
1433 }
1434
b8cd5251 1435 return reconstructor;
c757bafd 1436}
1437
2257f27e 1438//_____________________________________________________________________________
1439Bool_t AliReconstruction::CreateVertexer()
1440{
1441// create the vertexer
1442
b8cd5251 1443 fVertexer = NULL;
1444 AliReconstructor* itsReconstructor = GetReconstructor(0);
59697224 1445 if (itsReconstructor) {
b8cd5251 1446 fVertexer = itsReconstructor->CreateVertexer(fRunLoader);
2257f27e 1447 }
b8cd5251 1448 if (!fVertexer) {
815c2b38 1449 AliWarning("couldn't create a vertexer for ITS");
2257f27e 1450 if (fStopOnError) return kFALSE;
1451 }
1452
1453 return kTRUE;
1454}
1455
24f7a148 1456//_____________________________________________________________________________
b8cd5251 1457Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
24f7a148 1458{
f08fc9f5 1459// create the trackers
24f7a148 1460
b8cd5251 1461 TString detStr = detectors;
1462 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1463 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1464 AliReconstructor* reconstructor = GetReconstructor(iDet);
1465 if (!reconstructor) continue;
1466 TString detName = fgkDetectorName[iDet];
1f46a9ae 1467 if (detName == "HLT") {
1468 fRunHLTTracking = kTRUE;
1469 continue;
1470 }
f08fc9f5 1471
1472 fTracker[iDet] = reconstructor->CreateTracker(fRunLoader);
1473 if (!fTracker[iDet] && (iDet < 7)) {
1474 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
8250d5f5 1475 if (fStopOnError) return kFALSE;
1476 }
1477 }
1478
24f7a148 1479 return kTRUE;
1480}
1481
e583c30d 1482//_____________________________________________________________________________
b26c3770 1483void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
e583c30d 1484{
1485// delete trackers and the run loader and close and delete the file
1486
b8cd5251 1487 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1488 delete fReconstructor[iDet];
1489 fReconstructor[iDet] = NULL;
1490 fLoader[iDet] = NULL;
1491 delete fTracker[iDet];
1492 fTracker[iDet] = NULL;
1493 }
1494 delete fVertexer;
1495 fVertexer = NULL;
e583c30d 1496
1497 delete fRunLoader;
1498 fRunLoader = NULL;
b649205a 1499 delete fRawReader;
1500 fRawReader = NULL;
e583c30d 1501
1502 if (file) {
1503 file->Close();
1504 delete file;
1505 }
b26c3770 1506
1507 if (fileOld) {
1508 fileOld->Close();
1509 delete fileOld;
1510 gSystem->Unlink("AliESDs.old.root");
1511 }
e583c30d 1512}
24f7a148 1513
1514
1515//_____________________________________________________________________________
1516Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
1517{
1518// read the ESD event from a file
1519
1520 if (!esd) return kFALSE;
1521 char fileName[256];
1522 sprintf(fileName, "ESD_%d.%d_%s.root",
1523 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
1524 if (gSystem->AccessPathName(fileName)) return kFALSE;
1525
f3a97c86 1526 AliInfo(Form("reading ESD from file %s", fileName));
815c2b38 1527 AliDebug(1, Form("reading ESD from file %s", fileName));
24f7a148 1528 TFile* file = TFile::Open(fileName);
1529 if (!file || !file->IsOpen()) {
815c2b38 1530 AliError(Form("opening %s failed", fileName));
24f7a148 1531 delete file;
1532 return kFALSE;
1533 }
1534
1535 gROOT->cd();
1536 delete esd;
1537 esd = (AliESD*) file->Get("ESD");
1538 file->Close();
1539 delete file;
1540 return kTRUE;
1541}
1542
1543//_____________________________________________________________________________
1544void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
1545{
1546// write the ESD event to a file
1547
1548 if (!esd) return;
1549 char fileName[256];
1550 sprintf(fileName, "ESD_%d.%d_%s.root",
1551 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
1552
815c2b38 1553 AliDebug(1, Form("writing ESD to file %s", fileName));
24f7a148 1554 TFile* file = TFile::Open(fileName, "recreate");
1555 if (!file || !file->IsOpen()) {
815c2b38 1556 AliError(Form("opening %s failed", fileName));
24f7a148 1557 } else {
1558 esd->Write("ESD");
1559 file->Close();
1560 }
1561 delete file;
1562}
f3a97c86 1563
1564
1565
1566
1567//_____________________________________________________________________________
1568void AliReconstruction::CreateTag(TFile* file)
1569{
2bdb9d38 1570 /////////////
1571 //muon code//
1572 ////////////
56982dd1 1573 Double_t fMUONMASS = 0.105658369;
2bdb9d38 1574 //Variables
56982dd1 1575 Double_t fX,fY,fZ ;
1576 Double_t fThetaX, fThetaY, fPyz, fChisquare;
1577 Double_t fPxRec,fPyRec, fPzRec, fEnergy;
1578 Int_t fCharge;
1579 TLorentzVector fEPvector;
1580
1581 Float_t fZVertexCut = 10.0;
1582 Float_t fRhoVertexCut = 2.0;
1583
1584 Float_t fLowPtCut = 1.0;
1585 Float_t fHighPtCut = 3.0;
1586 Float_t fVeryHighPtCut = 10.0;
2bdb9d38 1587 ////////////
1588
1589 Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
1590
089bf903 1591 // Creates the tags for all the events in a given ESD file
4302e20f 1592 Int_t ntrack;
089bf903 1593 Int_t nProtons, nKaons, nPions, nMuons, nElectrons;
1594 Int_t nPos, nNeg, nNeutr;
1595 Int_t nK0s, nNeutrons, nPi0s, nGamas;
1596 Int_t nCh1GeV, nCh3GeV, nCh10GeV;
1597 Int_t nMu1GeV, nMu3GeV, nMu10GeV;
1598 Int_t nEl1GeV, nEl3GeV, nEl10GeV;
1599 Float_t maxPt = .0, meanPt = .0, totalP = .0;
d1a50cb5 1600 Int_t fVertexflag;
1387d165 1601 Int_t iRunNumber = 0;
547d75a6 1602 TString fVertexName("default");
4302e20f 1603
f3a97c86 1604 AliRunTag *tag = new AliRunTag();
f3a97c86 1605 AliEventTag *evTag = new AliEventTag();
1606 TTree ttag("T","A Tree with event tags");
2bdb9d38 1607 TBranch * btag = ttag.Branch("AliTAG", &tag);
f3a97c86 1608 btag->SetCompressionLevel(9);
2bdb9d38 1609
f3a97c86 1610 AliInfo(Form("Creating the tags......."));
1611
1612 if (!file || !file->IsOpen()) {
1613 AliError(Form("opening failed"));
1614 delete file;
1615 return ;
2bdb9d38 1616 }
1617 Int_t firstEvent = 0,lastEvent = 0;
f3a97c86 1618 TTree *t = (TTree*) file->Get("esdTree");
1619 TBranch * b = t->GetBranch("ESD");
1620 AliESD *esd = 0;
1621 b->SetAddress(&esd);
1387d165 1622
1623 b->GetEntry(0);
1624 Int_t iInitRunNumber = esd->GetRunNumber();
2bdb9d38 1625
089bf903 1626 Int_t iNumberOfEvents = b->GetEntries();
2bdb9d38 1627 for (Int_t iEventNumber = 0; iEventNumber < iNumberOfEvents; iEventNumber++) {
1628 ntrack = 0;
1629 nPos = 0;
1630 nNeg = 0;
1631 nNeutr =0;
1632 nK0s = 0;
1633 nNeutrons = 0;
1634 nPi0s = 0;
1635 nGamas = 0;
1636 nProtons = 0;
1637 nKaons = 0;
1638 nPions = 0;
1639 nMuons = 0;
1640 nElectrons = 0;
1641 nCh1GeV = 0;
1642 nCh3GeV = 0;
1643 nCh10GeV = 0;
1644 nMu1GeV = 0;
1645 nMu3GeV = 0;
1646 nMu10GeV = 0;
1647 nEl1GeV = 0;
1648 nEl3GeV = 0;
1649 nEl10GeV = 0;
1650 maxPt = .0;
1651 meanPt = .0;
1652 totalP = .0;
547d75a6 1653 fVertexflag = 0;
d1a50cb5 1654
2bdb9d38 1655 b->GetEntry(iEventNumber);
1387d165 1656 iRunNumber = esd->GetRunNumber();
1657 if(iRunNumber != iInitRunNumber) AliFatal("Inconsistency of run numbers in the AliESD!!!");
2bdb9d38 1658 const AliESDVertex * vertexIn = esd->GetVertex();
547d75a6 1659 if (!vertexIn) AliError("ESD has not defined vertex.");
1660 if (vertexIn) fVertexName = vertexIn->GetName();
1661 if(fVertexName != "default") fVertexflag = 1;
2bdb9d38 1662 for (Int_t iTrackNumber = 0; iTrackNumber < esd->GetNumberOfTracks(); iTrackNumber++) {
1663 AliESDtrack * esdTrack = esd->GetTrack(iTrackNumber);
1664 UInt_t status = esdTrack->GetStatus();
f3a97c86 1665
2bdb9d38 1666 //select only tracks with ITS refit
1667 if ((status&AliESDtrack::kITSrefit)==0) continue;
1668 //select only tracks with TPC refit
1669 if ((status&AliESDtrack::kTPCrefit)==0) continue;
4302e20f 1670
2bdb9d38 1671 //select only tracks with the "combined PID"
1672 if ((status&AliESDtrack::kESDpid)==0) continue;
1673 Double_t p[3];
1674 esdTrack->GetPxPyPz(p);
1675 Double_t momentum = sqrt(pow(p[0],2) + pow(p[1],2) + pow(p[2],2));
1676 Double_t fPt = sqrt(pow(p[0],2) + pow(p[1],2));
1677 totalP += momentum;
1678 meanPt += fPt;
1679 if(fPt > maxPt) maxPt = fPt;
4302e20f 1680
2bdb9d38 1681 if(esdTrack->GetSign() > 0) {
1682 nPos++;
56982dd1 1683 if(fPt > fLowPtCut) nCh1GeV++;
1684 if(fPt > fHighPtCut) nCh3GeV++;
1685 if(fPt > fVeryHighPtCut) nCh10GeV++;
2bdb9d38 1686 }
1687 if(esdTrack->GetSign() < 0) {
1688 nNeg++;
56982dd1 1689 if(fPt > fLowPtCut) nCh1GeV++;
1690 if(fPt > fHighPtCut) nCh3GeV++;
1691 if(fPt > fVeryHighPtCut) nCh10GeV++;
2bdb9d38 1692 }
1693 if(esdTrack->GetSign() == 0) nNeutr++;
4302e20f 1694
2bdb9d38 1695 //PID
1696 Double_t prob[5];
1697 esdTrack->GetESDpid(prob);
4302e20f 1698
2bdb9d38 1699 Double_t rcc = 0.0;
1700 for(Int_t i = 0; i < AliPID::kSPECIES; i++) rcc += prob[i]*partFrac[i];
1701 if(rcc == 0.0) continue;
1702 //Bayes' formula
1703 Double_t w[5];
1704 for(Int_t i = 0; i < AliPID::kSPECIES; i++) w[i] = prob[i]*partFrac[i]/rcc;
4302e20f 1705
2bdb9d38 1706 //protons
1707 if ((w[4]>w[3])&&(w[4]>w[2])&&(w[4]>w[1])&&(w[4]>w[0])) nProtons++;
1708 //kaons
1709 if ((w[3]>w[4])&&(w[3]>w[2])&&(w[3]>w[1])&&(w[3]>w[0])) nKaons++;
1710 //pions
1711 if ((w[2]>w[4])&&(w[2]>w[3])&&(w[2]>w[1])&&(w[2]>w[0])) nPions++;
1712 //electrons
1713 if ((w[0]>w[4])&&(w[0]>w[3])&&(w[0]>w[2])&&(w[0]>w[1])) {
1714 nElectrons++;
56982dd1 1715 if(fPt > fLowPtCut) nEl1GeV++;
1716 if(fPt > fHighPtCut) nEl3GeV++;
1717 if(fPt > fVeryHighPtCut) nEl10GeV++;
2bdb9d38 1718 }
1719 ntrack++;
1720 }//track loop
1721
1722 /////////////
1723 //muon code//
1724 ////////////
1725 Int_t nMuonTracks = esd->GetNumberOfMuonTracks();
1726 // loop over all reconstructed tracks (also first track of combination)
1727 for (Int_t iTrack = 0; iTrack < nMuonTracks; iTrack++) {
1728 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
1729 if (muonTrack == 0x0) continue;
4302e20f 1730
2bdb9d38 1731 // Coordinates at vertex
56982dd1 1732 fZ = muonTrack->GetZ();
1733 fY = muonTrack->GetBendingCoor();
1734 fX = muonTrack->GetNonBendingCoor();
4302e20f 1735
56982dd1 1736 fThetaX = muonTrack->GetThetaX();
1737 fThetaY = muonTrack->GetThetaY();
4302e20f 1738
56982dd1 1739 fPyz = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
1740 fPzRec = - fPyz / TMath::Sqrt(1.0 + TMath::Tan(fThetaY)*TMath::Tan(fThetaY));
1741 fPxRec = fPzRec * TMath::Tan(fThetaX);
1742 fPyRec = fPzRec * TMath::Tan(fThetaY);
1743 fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
f3a97c86 1744
2bdb9d38 1745 //ChiSquare of the track if needed
56982dd1 1746 fChisquare = muonTrack->GetChi2()/(2.0 * muonTrack->GetNHit() - 5);
1747 fEnergy = TMath::Sqrt(fMUONMASS * fMUONMASS + fPxRec * fPxRec + fPyRec * fPyRec + fPzRec * fPzRec);
1748 fEPvector.SetPxPyPzE(fPxRec, fPyRec, fPzRec, fEnergy);
4302e20f 1749
2bdb9d38 1750 // total number of muons inside a vertex cut
56982dd1 1751 if((TMath::Abs(fZ)<fZVertexCut) && (TMath::Sqrt(fY*fY+fX*fX)<fRhoVertexCut)) {
2bdb9d38 1752 nMuons++;
56982dd1 1753 if(fEPvector.Pt() > fLowPtCut) {
2bdb9d38 1754 nMu1GeV++;
56982dd1 1755 if(fEPvector.Pt() > fHighPtCut) {
2bdb9d38 1756 nMu3GeV++;
56982dd1 1757 if (fEPvector.Pt() > fVeryHighPtCut) {
2bdb9d38 1758 nMu10GeV++;
1759 }
1760 }
1761 }
1762 }
1763 }//muon track loop
1764
1765 // Fill the event tags
1766 if(ntrack != 0)
1767 meanPt = meanPt/ntrack;
1768
1769 evTag->SetEventId(iEventNumber+1);
547d75a6 1770 if (vertexIn) {
1771 evTag->SetVertexX(vertexIn->GetXv());
1772 evTag->SetVertexY(vertexIn->GetYv());
1773 evTag->SetVertexZ(vertexIn->GetZv());
1774 evTag->SetVertexZError(vertexIn->GetZRes());
1775 }
d1a50cb5 1776 evTag->SetVertexFlag(fVertexflag);
1777
2bdb9d38 1778 evTag->SetT0VertexZ(esd->GetT0zVertex());
1779
8bd8ac26 1780 evTag->SetTriggerMask(esd->GetTriggerMask());
1781 evTag->SetTriggerCluster(esd->GetTriggerCluster());
2bdb9d38 1782
32a5cab4 1783 evTag->SetZDCNeutron1Energy(esd->GetZDCN1Energy());
1784 evTag->SetZDCProton1Energy(esd->GetZDCP1Energy());
1785 evTag->SetZDCNeutron2Energy(esd->GetZDCN2Energy());
1786 evTag->SetZDCProton2Energy(esd->GetZDCP2Energy());
2bdb9d38 1787 evTag->SetZDCEMEnergy(esd->GetZDCEMEnergy());
1788 evTag->SetNumOfParticipants(esd->GetZDCParticipants());
1789
1790
1791 evTag->SetNumOfTracks(esd->GetNumberOfTracks());
1792 evTag->SetNumOfPosTracks(nPos);
1793 evTag->SetNumOfNegTracks(nNeg);
1794 evTag->SetNumOfNeutrTracks(nNeutr);
1795
1796 evTag->SetNumOfV0s(esd->GetNumberOfV0s());
1797 evTag->SetNumOfCascades(esd->GetNumberOfCascades());
1798 evTag->SetNumOfKinks(esd->GetNumberOfKinks());
1799 evTag->SetNumOfPMDTracks(esd->GetNumberOfPmdTracks());
1800
1801 evTag->SetNumOfProtons(nProtons);
1802 evTag->SetNumOfKaons(nKaons);
1803 evTag->SetNumOfPions(nPions);
1804 evTag->SetNumOfMuons(nMuons);
1805 evTag->SetNumOfElectrons(nElectrons);
1806 evTag->SetNumOfPhotons(nGamas);
1807 evTag->SetNumOfPi0s(nPi0s);
1808 evTag->SetNumOfNeutrons(nNeutrons);
1809 evTag->SetNumOfKaon0s(nK0s);
1810
1811 evTag->SetNumOfChargedAbove1GeV(nCh1GeV);
1812 evTag->SetNumOfChargedAbove3GeV(nCh3GeV);
1813 evTag->SetNumOfChargedAbove10GeV(nCh10GeV);
1814 evTag->SetNumOfMuonsAbove1GeV(nMu1GeV);
1815 evTag->SetNumOfMuonsAbove3GeV(nMu3GeV);
1816 evTag->SetNumOfMuonsAbove10GeV(nMu10GeV);
1817 evTag->SetNumOfElectronsAbove1GeV(nEl1GeV);
1818 evTag->SetNumOfElectronsAbove3GeV(nEl3GeV);
1819 evTag->SetNumOfElectronsAbove10GeV(nEl10GeV);
1820
85c60a8e 1821 evTag->SetNumOfPHOSClusters(esd->GetNumberOfPHOSClusters());
1822 evTag->SetNumOfEMCALClusters(esd->GetNumberOfEMCALClusters());
2bdb9d38 1823
1824 evTag->SetTotalMomentum(totalP);
1825 evTag->SetMeanPt(meanPt);
1826 evTag->SetMaxPt(maxPt);
1827
1387d165 1828 tag->SetRunId(iInitRunNumber);
2bdb9d38 1829 tag->AddEventTag(*evTag);
1830 }
089bf903 1831 lastEvent = iNumberOfEvents;
f3a97c86 1832
1833 ttag.Fill();
1834 tag->Clear();
1835
1836 char fileName[256];
1837 sprintf(fileName, "Run%d.Event%d_%d.ESD.tag.root",
1838 tag->GetRunId(),firstEvent,lastEvent );
1839 AliInfo(Form("writing tags to file %s", fileName));
1840 AliDebug(1, Form("writing tags to file %s", fileName));
1841
1842 TFile* ftag = TFile::Open(fileName, "recreate");
1843 ftag->cd();
1844 ttag.Write();
1845 ftag->Close();
1846 file->cd();
1847 delete tag;
f3a97c86 1848 delete evTag;
1849}
1850
98937d93 1851void AliReconstruction::WriteAlignmentData(AliESD* esd)
1852{
1853 // Write space-points which are then used in the alignment procedures
1854 // For the moment only ITS, TRD and TPC
1855
1856 // Load TOF clusters
d528ee75 1857 if (fTracker[3]){
1858 fLoader[3]->LoadRecPoints("read");
1859 TTree* tree = fLoader[3]->TreeR();
1860 if (!tree) {
1861 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
1862 return;
1863 }
1864 fTracker[3]->LoadClusters(tree);
98937d93 1865 }
98937d93 1866 Int_t ntracks = esd->GetNumberOfTracks();
1867 for (Int_t itrack = 0; itrack < ntracks; itrack++)
1868 {
1869 AliESDtrack *track = esd->GetTrack(itrack);
1870 Int_t nsp = 0;
ef7253ac 1871 Int_t idx[200];
98937d93 1872 for (Int_t iDet = 3; iDet >= 0; iDet--)
1873 nsp += track->GetNcls(iDet);
1874 if (nsp) {
1875 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
1876 track->SetTrackPointArray(sp);
1877 Int_t isptrack = 0;
1878 for (Int_t iDet = 3; iDet >= 0; iDet--) {
1879 AliTracker *tracker = fTracker[iDet];
1880 if (!tracker) continue;
1881 Int_t nspdet = track->GetNcls(iDet);
98937d93 1882 if (nspdet <= 0) continue;
1883 track->GetClusters(iDet,idx);
1884 AliTrackPoint p;
1885 Int_t isp = 0;
1886 Int_t isp2 = 0;
1887 while (isp < nspdet) {
1888 Bool_t isvalid = tracker->GetTrackPoint(idx[isp2],p); isp2++;
160db090 1889 const Int_t kNTPCmax = 159;
1890 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
98937d93 1891 if (!isvalid) continue;
1892 sp->AddPoint(isptrack,&p); isptrack++; isp++;
1893 }
98937d93 1894 }
1895 }
1896 }
d528ee75 1897 if (fTracker[3]){
1898 fTracker[3]->UnloadClusters();
1899 fLoader[3]->UnloadRecPoints();
1900 }
98937d93 1901}