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