]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFSDigits2Digits.C
Moving to the new VMC naming convention
[u/mrichter/AliRoot.git] / TOF / AliTOFSDigits2Digits.C
CommitLineData
c630aafd 1Int_t AliTOFSDigits2Digits(Int_t nev=5) {
ca4508a6 2
c630aafd 3 // Adapted to the NewIO by I.Belikov (Jouri.Belikov@cern.ch)
ca4508a6 4
5 // number 3 is a legacy from AliDigit object
ca4508a6 6
c630aafd 7 if (gAlice) {
8 delete gAlice->GetRunLoader();
9 delete gAlice;//if everything was OK here it is already NULL
10 gAlice = 0x0;
11 }
12
13 AliRunLoader* rl = AliRunLoader::Open("galice.root");
14 if (rl == 0x0) {
15 cerr<<"Can not open session"<<endl;
16 return 1;
17 }
18
19 if (rl->LoadgAlice()) {
20 cerr<<"Error occured while loading gAlice"<<endl;
21 return 1;
22 }
23
24 AliLoader *tofl = rl->GetLoader("TOFLoader");
25 if (tofl == 0x0) {
26 cerr<<"Can not get the TOF Loader"<<endl;
27 return 1;
28 }
29
30 gAlice=rl->GetAliRun();
31 if (!gAlice) {
32 cerr<<"Can't get gAlice !\n";
33 return 1;
34 }
35
36 tofl->LoadSDigits("read");
37 tofl->LoadDigits("recreate");
38
39 Int_t totndig=0; // total number of digits
40 Int_t tottracks=0; // total number of tracks contributing to totndig
41
42 if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
43
44 TClonesArray *fSDigits=new TClonesArray("AliTOFSDigit", 1000);
45 TClonesArray *fDigits =new TClonesArray("AliTOFdigit", 1000);
46 TClonesArray &da=*fDigits;
47
48 for (Int_t ievent = 0; ievent < nev; ievent++) {
49 rl->GetEvent(ievent);
50
51 TTree *sTree=tofl->TreeS();
52 if (sTree == 0) {
53 cerr<<"Can't get the sdigit tree !\n";
54 return 1;
55 }
56 TBranch *branch=sTree->GetBranch("TOF");
57 if (!branch) {
58 cerr<<"Cant' get the branch !\n";
59 return 1;
60 }
61 branch->SetAddress(&fSDigits);
ca4508a6 62
c630aafd 63 TTree *dTree=tofl->TreeD();
64 if (dTree == 0) {
65 tofl->MakeTree("D");
66 dTree=tofl->TreeD();
67 }
68 branch=dTree->GetBranch("TOF");
69 if (!branch) dTree->Branch("TOF",&fDigits);
70 else branch->SetAddress(&fDigits);
71
72 Int_t nEntries = sTree->GetEntries();
ca4508a6 73 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
c630aafd 74 sTree->GetEvent(iEntry);
ca4508a6 75
ca4508a6 76 Int_t ndig = fSDigits->GetEntriesFast();
c630aafd 77 cout << "----------------<AliTOFSDigits2Digits>----------------" << endl;
ca4508a6 78 cout << "Found " << ndig << " TOF SDigits for event " << ievent << endl;
c630aafd 79 cout << "------------------------------------------------------" << endl;
ca4508a6 80
ca4508a6 81 for (Int_t k = 0; k < ndig; k++) {
ca4508a6 82 Int_t vol[5]; // location for a digit
83
84 // Get the information for this digit
c630aafd 85 AliTOFSDigit *tofsdigit = (AliTOFSDigit *)fSDigits->UncheckedAt(k);
ca4508a6 86
87 Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
88 // for current sdigit
89
90 // TOF sdigit volumes (always the same for all slots)
91 Int_t sector = tofsdigit->GetSector(); // range [1-18]
92 Int_t plate = tofsdigit->GetPlate(); // range [1- 5]
93 Int_t strip = tofsdigit->GetStrip(); // range [1-20]
94 Int_t padz = tofsdigit->GetPadz(); // range [1- 2]
95 Int_t padx = tofsdigit->GetPadx(); // range [1-48]
96
97 vol[0] = sector;
98 vol[1] = plate;
99 vol[2] = strip;
100 vol[3] = padx;
101 vol[4] = padz;
102
103 //--------------------- QA section ----------------------
104 // in the while, I perform QA
c630aafd 105 Bool_t isSDigitBad = (sector<1 || sector>18 ||
106 plate<1 || plate >5 ||
107 padz<1 || padz>2 ||
108 padx<1 || padx>48);
ca4508a6 109
110 if (isSDigitBad) {
111 cout << "<AliTOFSDigits2Digits> strange sdigit found" << endl;
c630aafd 112 return 3;
ca4508a6 113 }
114 //-------------------------------------------------------
115
ca4508a6 116 // start loop on number of slots for current sdigit
117 for (Int_t islot = 0; islot < nslot; islot++) {
118 Float_t digit[2]; // TOF digit variables
c630aafd 119 const Int_t kMAXDIGITS = 3;
120 Int_t tracknum[kMAXDIGITS];//contributing tracks for the current slot
ca4508a6 121
122 Float_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
123 Float_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
124
125 tracknum[0]=tofsdigit->GetTrack(islot,0);
126 tracknum[1]=tofsdigit->GetTrack(islot,1);
127 tracknum[2]=tofsdigit->GetTrack(islot,2);
128
129 for (Int_t i = 0; i < kMAXDIGITS; i++) {
130 tottracks++;
131 // search for the first empty location
132 if(tracknum[i]==-1){
133 tottracks--;
134 break;
135 }
136 }
137
138 // adding a TOF digit for each slot
c630aafd 139 {
140 Int_t ndigits=da.GetEntriesFast();
141 //cerr<<ndigits<<" ndig \n";
142 new (da[ndigits]) AliTOFdigit(tracknum, vol, digit);
143 }
ca4508a6 144 totndig++;
145 }
146
147
148 } // end loop on sdigits
c630aafd 149 fSDigits->Clear();
ca4508a6 150
151 } // end loop on entries
152
c630aafd 153 dTree->Fill();
154 tofl->WriteDigits("OVERWRITE");
ca4508a6 155
c630aafd 156 // free used memory
157 fDigits->Clear();
ca4508a6 158
159 } // end loop on events
160
c630aafd 161 delete fSDigits;
162 delete fDigits;
163
ca4508a6 164 cout << "----------------------------------------------------------" << endl;
165 cout << "<AliTOFSDigits2Digits> Summary" << endl;
166 cout << "contributing tracks to " << totndig << " digits: " << tottracks << endl;
167 cout << "----------------------------------------------------------" << endl;
168
c630aafd 169 return 0;
ca4508a6 170
171
172}