1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // A. De Falco, H. Woehri, INFN Cagliari, July 2006
19 // This macro reads the generation/reconstruction files in the
20 // input directory (dirname),
21 // builds a tree that contains, for each event, an array of muons and dimuons
22 // (TClonesArrays of AliMUONTrackLight and AliMUONPairLight objects)
23 // and writes the tree in an output file (outFileName)
24 // Note that if the path for the output file is not explicitly specified,
25 // it will be written in the directory containing the generation/reconstruction
28 #include <TClonesArray.h>
32 #include <TLorentzVector.h>
33 #include <TParticle.h>
35 #include <AliMUONRecoCheck.h>
36 #include <AliMUONTrack.h>
37 #include <AliMUONTrackParam.h>
38 #include <AliRunLoader.h>
39 #include "AliMUONTrackLight.h"
40 #include "AliMUONPairLight.h"
41 #include "AliMUONTrackExtrap.h"
43 void DecodeRecoCocktail(char* dirname=".", char* outFileName = "MuonLight.root"){
44 const char *startingDir = gSystem->pwd();
46 TFile *fout = new TFile(outFileName,"recreate");
48 AliMUONRecoCheck *rc = new AliMUONRecoCheck("galice.root");
49 AliRunLoader *runLoader = rc->GetRunLoader();
51 TClonesArray *muonArray = new TClonesArray("AliMUONTrackLight",10);
52 TClonesArray *dimuonArray = new TClonesArray("AliMUONPairLight",10);
53 TTree *treeOut = new TTree("tree","tree");
56 treeOut->Branch("muons",&muonArray);
57 treeOut->Branch("dimuons",&dimuonArray);
59 runLoader->LoadKinematics("READ");;
60 Int_t nev = runLoader->GetNumberOfEvents();
63 for(Int_t ievent = 0; ievent < nev; ievent++){ // loop over events
64 printf ("Event %d of %d\n",ievent,nev);
65 muonArray->Clear(); // clean muon and dimuon arrays
67 runLoader->GetEvent(ievent);
69 rc->MakeTrackRef(); // make reconstructable tracks
71 TClonesArray * trackRecoArray = rc->GetTrackReco();
72 TClonesArray * trackRefArray = rc->GetMuonTrackRef();
73 Int_t nTrackReco = trackRecoArray->GetEntriesFast();
75 for (Int_t itrRec = 0; itrRec<nTrackReco; itrRec++) {
76 // assign parameters concerning the reconstructed tracks
77 AliMUONTrack *trackReco = (AliMUONTrack *)trackRecoArray->At(itrRec);
78 AliMUONTrackLight muLight;
79 AliMUONTrackParam trPar(*((AliMUONTrackParam*) (trackReco->GetTrackParamAtHit()->First())));
80 AliMUONTrackExtrap::ExtrapToVertex(&trPar,0.,0.,0.);
81 muLight.SetCharge(Int_t(TMath::Sign(1.,trPar.GetInverseBendingMomentum())));
82 muLight.SetPxPyPz(trPar.Px(),trPar.Py(), trPar.Pz());
83 muLight.SetTriggered(trackReco->GetMatchTrigger());
84 Double_t xyz[3] = { trPar.GetNonBendingCoor(),
85 trPar.GetBendingCoor(),
87 muLight.SetVertex(xyz);
88 // find the reference track and store further information
89 TParticle *part = muLight.FindRefTrack(trackReco,trackRefArray,runLoader);
91 v.SetPxPyPzE(part->Px(), part->Py(), part->Pz(), part->Energy());
93 muLight.FillMuonHistory(runLoader, part);
94 muLight.PrintInfo("A");
95 //store the referenced track in the muonArray:
96 TClonesArray &muons = *muonArray;
97 new (muons[nreftracks++]) AliMUONTrackLight(muLight);
101 // now loop over muon pairs to build dimuons
102 Int_t nmuons = muonArray->GetEntriesFast();
104 for(Int_t itrRec1 = 0; itrRec1 < nmuons-1; itrRec1++) {
105 AliMUONTrackLight* mu1 = (AliMUONTrackLight*) muonArray->At(itrRec1);
106 for(Int_t itrRec2 = itrRec1+1; itrRec2 < nmuons; itrRec2++){
107 AliMUONTrackLight* mu2 = (AliMUONTrackLight*) muonArray->At(itrRec2);
108 AliMUONPairLight dimuLight;
109 dimuLight.SetMuons(*mu1, *mu2);
110 // if(dimuLight.GetCreationProcess() == 2){
111 // printf ("#dimuon = %d (%d, %d) \n", ndimuons, itrRec1, itrRec2);
112 // dimuLight.PrintInfo("A");
114 //store the referenced track in the dimuonArray:
115 TClonesArray &dimuons = *dimuonArray;
116 new (dimuons[ndimuons++]) AliMUONPairLight(dimuLight);
119 //Int_t ndimu2 = dimuonArray->GetEntriesFast();
120 // printf ("dimuonArray has %d entries\n",ndimu2);
121 // dimuonArray->Dump();
122 // printf ("filling tree\n");
125 // printf ("done\n");
130 gSystem->cd(startingDir);