7587f5a5 |
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 | |
b2a60966 |
16 | /* $Id$ */ |
5f20d3fb |
17 | |
7587f5a5 |
18 | //_________________________________________________________________________ |
5f20d3fb |
19 | // Implementation version v1 of PHOS Manager class |
a3dfe79c |
20 | //--- |
21 | // Layout EMC + PPSD has name GPS2: |
5f20d3fb |
22 | // Produces cumulated hits (no hits) and digits |
a3dfe79c |
23 | //--- |
24 | // Layout EMC + CPV has name IHEP: |
25 | // Produces hits for CPV, cumulated hits and digits |
5f20d3fb |
26 | //*-- Author: Yves Schutz (SUBATECH) |
b2a60966 |
27 | |
7587f5a5 |
28 | |
29 | // --- ROOT system --- |
bea63bea |
30 | |
31 | #include "TBRIK.h" |
32 | #include "TNode.h" |
7587f5a5 |
33 | #include "TRandom.h" |
94de3818 |
34 | #include "TTree.h" |
7587f5a5 |
35 | |
5f20d3fb |
36 | |
7587f5a5 |
37 | // --- Standard library --- |
38 | |
de9ec31b |
39 | #include <stdio.h> |
40 | #include <string.h> |
41 | #include <stdlib.h> |
42 | #include <strstream.h> |
7587f5a5 |
43 | |
44 | // --- AliRoot header files --- |
45 | |
46 | #include "AliPHOSv1.h" |
47 | #include "AliPHOSHit.h" |
48 | #include "AliPHOSDigit.h" |
bea63bea |
49 | #include "AliPHOSReconstructioner.h" |
7587f5a5 |
50 | #include "AliRun.h" |
51 | #include "AliConst.h" |
94de3818 |
52 | #include "AliMC.h" |
7587f5a5 |
53 | |
54 | ClassImp(AliPHOSv1) |
55 | |
bea63bea |
56 | //____________________________________________________________________________ |
57 | AliPHOSv1::AliPHOSv1() |
58 | { |
5f20d3fb |
59 | // ctor |
bea63bea |
60 | fNTmpHits = 0 ; |
fa412d9b |
61 | fTmpHits = 0 ; |
fc879520 |
62 | |
3d402178 |
63 | // Create an empty array of AliPHOSCPVModule to satisfy |
fa412d9b |
64 | // AliPHOSv1::Streamer when reading root file |
65 | |
3d402178 |
66 | if ( NULL==(fCPVModules=new TClonesArray("AliPHOSCPVModule",0)) ) { |
fa412d9b |
67 | Error("AliPHOSv1","Can not create array of CPV modules"); |
68 | exit(1); |
69 | } |
fc879520 |
70 | |
bea63bea |
71 | } |
72 | |
7587f5a5 |
73 | //____________________________________________________________________________ |
74 | AliPHOSv1::AliPHOSv1(const char *name, const char *title): |
e04976bd |
75 | AliPHOSv0(name,title) |
7587f5a5 |
76 | { |
5f20d3fb |
77 | // ctor : title is used to identify the layout |
fa412d9b |
78 | // GPS2 = 5 modules (EMC + PPSD) |
79 | // IHEP = 5 modules (EMC + CPV ) |
5f20d3fb |
80 | // We use 2 arrays of hits : |
81 | // |
82 | // - fHits (the "normal" one), which retains the hits associated with |
83 | // the current primary particle being tracked |
84 | // (this array is reset after each primary has been tracked). |
85 | // |
86 | // - fTmpHits, which retains all the hits of the current event. It |
87 | // is used for the digitization part. |
fa412d9b |
88 | |
5f20d3fb |
89 | fPinElectronicNoise = 0.010 ; |
b95514c0 |
90 | fDigitThreshold = 0.1 ; // 1 GeV |
5f20d3fb |
91 | |
92 | // We do not want to save in TreeH the raw hits |
93 | // But save the cumulated hits instead (need to create the branch myself) |
94 | // It is put in the Digit Tree because the TreeH is filled after each primary |
fa412d9b |
95 | // and the TreeD at the end of the event (branch is set in FinishEvent() ). |
5f20d3fb |
96 | |
97 | fTmpHits= new TClonesArray("AliPHOSHit",1000) ; |
98 | |
99 | fNTmpHits = fNhits = 0 ; |
100 | |
101 | fDigits = new TClonesArray("AliPHOSDigit",1000) ; |
102 | |
5f20d3fb |
103 | fIshunt = 1 ; // All hits are associated with primary particles |
fa412d9b |
104 | |
105 | // Create array of CPV modules for the IHEP's version of CPV |
106 | |
107 | if ( strcmp(fGeom->GetName(),"IHEP") == 0 ) { |
3d402178 |
108 | // Create array of AliPHOSCPVmodule of the size of PHOS modules number |
fa412d9b |
109 | |
3d402178 |
110 | if ( NULL==(fCPVModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNModules())) ) { |
fa412d9b |
111 | Error("AliPHOSv1","Can not create array of CPV modules"); |
112 | exit(1); |
113 | } |
114 | TClonesArray &lcpvmodule = *fCPVModules; |
3d402178 |
115 | for (Int_t i=0; i<fGeom->GetNModules(); i++) new(lcpvmodule[i]) AliPHOSCPVModule(); |
fa412d9b |
116 | } |
117 | else { |
3d402178 |
118 | // Create an empty array of AliPHOSCPVModule to satisfy |
fa412d9b |
119 | // AliPHOSv1::Streamer when writing root file |
120 | |
3d402178 |
121 | fCPVModules=new TClonesArray("AliPHOSCPVModule",0); |
fa412d9b |
122 | |
123 | } |
5f20d3fb |
124 | } |
125 | |
126 | //____________________________________________________________________________ |
127 | AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title): |
128 | AliPHOSv0(name,title) |
129 | { |
130 | // ctor : title is used to identify the layout |
131 | // GPS2 = 5 modules (EMC + PPSD) |
132 | // We use 2 arrays of hits : |
133 | // |
134 | // - fHits (the "normal" one), which retains the hits associated with |
135 | // the current primary particle being tracked |
136 | // (this array is reset after each primary has been tracked). |
137 | // |
138 | // - fTmpHits, which retains all the hits of the current event. It |
139 | // is used for the digitization part. |
140 | |
141 | fPinElectronicNoise = 0.010 ; |
142 | |
143 | // We do not want to save in TreeH the raw hits |
144 | //fHits = new TClonesArray("AliPHOSHit",100) ; |
145 | |
146 | fDigits = new TClonesArray("AliPHOSDigit",1000) ; |
147 | fTmpHits= new TClonesArray("AliPHOSHit",1000) ; |
148 | |
149 | fNTmpHits = fNhits = 0 ; |
150 | |
151 | fIshunt = 1 ; // All hits are associated with primary particles |
152 | |
153 | // gets an instance of the geometry parameters class |
154 | fGeom = AliPHOSGeometry::GetInstance(title, "") ; |
155 | |
156 | if (fGeom->IsInitialized() ) |
88bdfa12 |
157 | cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ; |
5f20d3fb |
158 | else |
fc879520 |
159 | cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ; |
5f20d3fb |
160 | |
161 | // Defining the PHOS Reconstructioner |
162 | |
163 | fReconstructioner = Reconstructioner ; |
164 | |
7587f5a5 |
165 | } |
b2a60966 |
166 | |
7587f5a5 |
167 | //____________________________________________________________________________ |
bea63bea |
168 | AliPHOSv1::~AliPHOSv1() |
b2a60966 |
169 | { |
bea63bea |
170 | // dtor |
5f20d3fb |
171 | |
8dfa469d |
172 | if ( fTmpHits) { |
173 | fTmpHits->Delete() ; |
174 | delete fTmpHits ; |
175 | fTmpHits = 0 ; |
176 | } |
5f20d3fb |
177 | |
fc879520 |
178 | // if ( fEmcRecPoints ) { |
179 | // fEmcRecPoints->Delete() ; |
180 | // delete fEmcRecPoints ; |
181 | // fEmcRecPoints = 0 ; |
182 | // } |
183 | |
184 | // if ( fPpsdRecPoints ) { |
185 | // fPpsdRecPoints->Delete() ; |
186 | // delete fPpsdRecPoints ; |
187 | // fPpsdRecPoints = 0 ; |
188 | // } |
5f20d3fb |
189 | |
fc879520 |
190 | // if ( fTrackSegments ) { |
191 | // fTrackSegments->Delete() ; |
192 | // delete fTrackSegments ; |
193 | // fTrackSegments = 0 ; |
194 | // } |
5f20d3fb |
195 | |
7587f5a5 |
196 | } |
197 | |
7587f5a5 |
198 | //____________________________________________________________________________ |
7eb9d12d |
199 | void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits, Int_t trackpid) |
bea63bea |
200 | { |
201 | // Add a hit to the hit list. |
5f20d3fb |
202 | // A PHOS hit is the sum of all hits in a single crystal |
203 | // or in a single PPSD gas cell |
bea63bea |
204 | |
5f20d3fb |
205 | Int_t hitCounter ; |
206 | TClonesArray <mphits = *fTmpHits ; |
bea63bea |
207 | AliPHOSHit *newHit ; |
5f20d3fb |
208 | AliPHOSHit *curHit ; |
209 | Bool_t deja = kFALSE ; |
bea63bea |
210 | |
5f20d3fb |
211 | // In any case, fills the fTmpHit TClonesArray (with "accumulated hits") |
bea63bea |
212 | |
7eb9d12d |
213 | newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid) ; |
bea63bea |
214 | |
5f20d3fb |
215 | // We do not want to save in TreeH the raw hits |
bea63bea |
216 | // TClonesArray &lhits = *fHits; |
bea63bea |
217 | |
5f20d3fb |
218 | for ( hitCounter = 0 ; hitCounter < fNTmpHits && !deja ; hitCounter++ ) { |
219 | curHit = (AliPHOSHit*) ltmphits[hitCounter] ; |
220 | if( *curHit == *newHit ) { |
221 | *curHit = *curHit + *newHit ; |
222 | deja = kTRUE ; |
223 | } |
224 | } |
225 | |
226 | if ( !deja ) { |
227 | new(ltmphits[fNTmpHits]) AliPHOSHit(*newHit) ; |
228 | fNTmpHits++ ; |
229 | } |
230 | |
bea63bea |
231 | // We do not want to save in TreeH the raw hits |
232 | // new(lhits[fNhits]) AliPHOSHit(*newHit) ; |
233 | // fNhits++ ; |
234 | |
235 | // Please note that the fTmpHits array must survive up to the |
236 | // end of the events, so it does not appear e.g. in ResetHits() ( |
237 | // which is called at the end of each primary). |
238 | |
239 | delete newHit; |
240 | |
241 | } |
242 | |
5f20d3fb |
243 | //___________________________________________________________________________ |
244 | Int_t AliPHOSv1::Digitize(Float_t Energy) |
245 | { |
246 | // Applies the energy calibration |
247 | |
248 | Float_t fB = 100000000. ; |
249 | Float_t fA = 0. ; |
250 | Int_t chan = Int_t(fA + Energy*fB ) ; |
251 | return chan ; |
252 | } |
bea63bea |
253 | |
254 | //___________________________________________________________________________ |
255 | void AliPHOSv1::FinishEvent() |
256 | { |
257 | // Makes the digits from the sum of summed hit in a single crystal or PPSD gas cell |
258 | // Adds to the energy the electronic noise |
259 | // Keeps digits with energy above fDigitThreshold |
260 | |
261 | // Save the cumulated hits instead of raw hits (need to create the branch myself) |
262 | // It is put in the Digit Tree because the TreeH is filled after each primary |
263 | // and the TreeD at the end of the event. |
264 | |
5f20d3fb |
265 | |
bea63bea |
266 | Int_t i ; |
267 | Int_t relid[4]; |
268 | Int_t j ; |
269 | TClonesArray &lDigits = *fDigits ; |
fa412d9b |
270 | AliPHOSHit * hit ; |
bea63bea |
271 | AliPHOSDigit * newdigit ; |
272 | AliPHOSDigit * curdigit ; |
273 | Bool_t deja = kFALSE ; |
274 | |
5f20d3fb |
275 | for ( i = 0 ; i < fNTmpHits ; i++ ) { |
276 | hit = (AliPHOSHit*)fTmpHits->At(i) ; |
4a2ca5e9 |
277 | |
278 | // Assign primary number only if contribution is significant |
279 | if( hit->GetEnergy() > fDigitThreshold) |
280 | newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ; |
281 | else |
fa412d9b |
282 | newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ; |
bea63bea |
283 | deja =kFALSE ; |
284 | for ( j = 0 ; j < fNdigits ; j++) { |
285 | curdigit = (AliPHOSDigit*) lDigits[j] ; |
286 | if ( *curdigit == *newdigit) { |
287 | *curdigit = *curdigit + *newdigit ; |
288 | deja = kTRUE ; |
289 | } |
290 | } |
291 | if ( !deja ) { |
292 | new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ; |
293 | fNdigits++ ; |
294 | } |
295 | |
296 | delete newdigit ; |
297 | } |
298 | |
299 | // Noise induced by the PIN diode of the PbWO crystals |
300 | |
301 | Float_t energyandnoise ; |
302 | for ( i = 0 ; i < fNdigits ; i++ ) { |
303 | newdigit = (AliPHOSDigit * ) fDigits->At(i) ; |
fa412d9b |
304 | |
bea63bea |
305 | fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ; |
306 | |
307 | if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals) |
308 | energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ; |
309 | |
310 | if (energyandnoise < 0 ) |
311 | energyandnoise = 0 ; |
312 | |
313 | if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list |
314 | fDigits->RemoveAt(i) ; |
315 | } |
316 | } |
317 | |
318 | fDigits->Compress() ; |
319 | |
320 | fNdigits = fDigits->GetEntries() ; |
321 | for (i = 0 ; i < fNdigits ; i++) { |
322 | newdigit = (AliPHOSDigit *) fDigits->At(i) ; |
323 | newdigit->SetIndexInList(i) ; |
fa412d9b |
324 | |
325 | // fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ; |
326 | // printf("FinishEvent(): relid=(%d,%d,%d,%d) Amp=%d\n", |
327 | // relid[0],relid[1],relid[2],relid[3], newdigit->GetAmp()); |
bea63bea |
328 | } |
fa412d9b |
329 | |
5f20d3fb |
330 | } |
331 | |
332 | //___________________________________________________________________________ |
333 | void AliPHOSv1::MakeBranch(Option_t* opt) |
334 | { |
335 | // Create new branche in the current Root Tree in the digit Tree |
5f20d3fb |
336 | AliDetector::MakeBranch(opt) ; |
337 | |
338 | char branchname[10]; |
7eb9d12d |
339 | |
340 | // Create new branche PHOSCH in the current Root Tree in the Hit Tree for accumulated Hits |
341 | if ( ! (gAlice->IsLegoRun()) ) { // only when not in lego plot mode |
342 | char *cdH = strstr(opt,"H"); |
343 | if ( fTmpHits && gAlice->TreeH() && cdH) { |
344 | sprintf(branchname, "%sCH", GetName()) ; |
345 | gAlice->TreeH()->Branch(branchname, &fTmpHits, fBufferSize) ; |
346 | } |
347 | } |
5f20d3fb |
348 | sprintf(branchname,"%s",GetName()); |
349 | char *cdD = strstr(opt,"D"); |
350 | if (fDigits && gAlice->TreeD() && cdD) { |
351 | gAlice->TreeD()->Branch(branchname, &fDigits, fBufferSize); |
352 | } |
353 | |
bea63bea |
354 | |
fa412d9b |
355 | // Create new branches CPV<i> for hits in CPV modules for IHEP geometry |
356 | // Yuri Kharlov, 28 September 2000. |
357 | |
358 | if ( strcmp(fGeom->GetName(),"IHEP") == 0 ) { |
359 | for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetCPVModule(i).MakeBranch(i+1); |
360 | } |
361 | |
bea63bea |
362 | } |
363 | |
5f20d3fb |
364 | //_____________________________________________________________________________ |
365 | void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner) |
366 | { |
367 | // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and |
368 | // 2. Creates TreeR with a branch for each list |
369 | // 3. Steers the reconstruction processes |
370 | // 4. Saves the 3 lists in TreeR |
371 | // 5. Write the Tree to File |
372 | |
373 | fReconstructioner = Reconstructioner ; |
374 | |
375 | char branchname[10] ; |
376 | |
377 | // 1. |
378 | |
379 | // gAlice->MakeTree("R") ; |
380 | Int_t splitlevel = 0 ; |
381 | |
fc879520 |
382 | fEmcRecPoints->Delete() ; |
5f20d3fb |
383 | |
384 | if ( fEmcRecPoints && gAlice->TreeR() ) { |
385 | sprintf(branchname,"%sEmcRP",GetName()) ; |
5f20d3fb |
386 | gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ; |
387 | } |
388 | |
fa412d9b |
389 | fPpsdRecPoints->Delete() ; |
5f20d3fb |
390 | |
5f20d3fb |
391 | if ( fPpsdRecPoints && gAlice->TreeR() ) { |
392 | sprintf(branchname,"%sPpsdRP",GetName()) ; |
5f20d3fb |
393 | gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ; |
394 | } |
395 | |
fc879520 |
396 | fTrackSegments->Delete() ; |
5f20d3fb |
397 | |
5f20d3fb |
398 | if ( fTrackSegments && gAlice->TreeR() ) { |
399 | sprintf(branchname,"%sTS",GetName()) ; |
400 | gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ; |
401 | } |
402 | |
fa412d9b |
403 | fRecParticles->Delete() ; |
fc879520 |
404 | |
fa412d9b |
405 | if (strcmp(fGeom->GetName(),"GPS2") == 0) { |
406 | if ( fRecParticles && gAlice->TreeR() ) { |
407 | sprintf(branchname,"%sRP",GetName()) ; |
408 | gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ; |
409 | } |
5f20d3fb |
410 | } |
411 | |
412 | // 3. |
fa412d9b |
413 | if (strcmp(fGeom->GetName(),"GPS2") == 0) |
414 | fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles); |
415 | else if (strcmp(fGeom->GetName(),"IHEP") == 0) |
416 | fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints); |
5f20d3fb |
417 | |
418 | // 4. Expand or Shrink the arrays to the proper size |
419 | |
420 | Int_t size ; |
421 | |
422 | size = fEmcRecPoints->GetEntries() ; |
423 | fEmcRecPoints->Expand(size) ; |
424 | |
425 | size = fPpsdRecPoints->GetEntries() ; |
426 | fPpsdRecPoints->Expand(size) ; |
427 | |
428 | size = fTrackSegments->GetEntries() ; |
429 | fTrackSegments->Expand(size) ; |
430 | |
431 | size = fRecParticles->GetEntries() ; |
432 | fRecParticles->Expand(size) ; |
433 | |
434 | gAlice->TreeR()->Fill() ; |
5f20d3fb |
435 | // 5. |
436 | |
4a2ca5e9 |
437 | gAlice->TreeR()->Write(0,TObject::kOverwrite) ; |
5f20d3fb |
438 | |
439 | // Deleting reconstructed objects |
440 | ResetReconstruction(); |
441 | |
442 | } |
443 | |
fa412d9b |
444 | //____________________________________________________________________________ |
445 | void AliPHOSv1::ResetHits() |
c4f224e7 |
446 | { |
fa412d9b |
447 | // Reset hit tree for CPV in IHEP geometry |
448 | // Yuri Kharlov, 28 September 2000 |
449 | |
c4f224e7 |
450 | |
a16726d7 |
451 | if ( fTmpHits ) { |
452 | fTmpHits->Clear() ; |
453 | fNTmpHits = 0; |
454 | } |
c4f224e7 |
455 | |
fa412d9b |
456 | AliDetector::ResetHits(); |
457 | if ( strcmp(fGeom->GetName(),"IHEP") == 0 ) { |
3d402178 |
458 | for (Int_t i=0; i<fGeom->GetNModules(); i++) ((AliPHOSCPVModule*)(*fCPVModules)[i]) -> Clear(); |
fa412d9b |
459 | } |
460 | } |
5f20d3fb |
461 | //____________________________________________________________________________ |
462 | void AliPHOSv1::ResetDigits() |
463 | { |
464 | // May sound strange, but cumulative hits are store in digits Tree |
465 | AliDetector::ResetDigits(); |
466 | if( fTmpHits ) { |
467 | fTmpHits->Delete(); |
468 | fNTmpHits = 0 ; |
469 | } |
470 | } |
471 | //____________________________________________________________________________ |
472 | void AliPHOSv1::ResetReconstruction() |
473 | { |
474 | // Deleting reconstructed objects |
475 | |
476 | if ( fEmcRecPoints ) fEmcRecPoints->Delete(); |
477 | if ( fPpsdRecPoints ) fPpsdRecPoints->Delete(); |
478 | if ( fTrackSegments ) fTrackSegments->Delete(); |
479 | if ( fRecParticles ) fRecParticles->Delete(); |
480 | |
481 | } |
5f20d3fb |
482 | |
483 | //____________________________________________________________________________ |
484 | void AliPHOSv1::SetTreeAddress() |
485 | { |
486 | // TBranch *branch; |
487 | AliPHOS::SetTreeAddress(); |
488 | |
489 | // //Branch address for TreeR: RecPpsdRecPoint |
490 | // TTree *treeR = gAlice->TreeR(); |
491 | // if ( treeR && fPpsdRecPoints ) { |
492 | // branch = treeR->GetBranch("PHOSPpsdRP"); |
493 | // if (branch) branch->SetAddress(&fPpsdRecPoints) ; |
fa412d9b |
494 | // } |
495 | |
496 | // Set branch address for the Hits Tree for hits in CPV modules for IHEP geometry |
497 | // Yuri Kharlov, 28 September 2000. |
498 | |
499 | if ( strcmp(fGeom->GetName(),"IHEP") == 0 ) { |
500 | for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetCPVModule(i).SetTreeAddress(i+1); |
501 | } |
502 | |
5f20d3fb |
503 | } |
504 | |
505 | //____________________________________________________________________________ |
506 | |
7587f5a5 |
507 | void AliPHOSv1::StepManager(void) |
508 | { |
b2a60966 |
509 | // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell |
b2a60966 |
510 | |
7587f5a5 |
511 | Int_t relid[4] ; // (box, layer, row, column) indices |
fa412d9b |
512 | Int_t absid ; // absolute cell ID number |
513 | Float_t xyze[4] ; // position wrt MRS and energy deposited |
514 | TLorentzVector pos ; // Lorentz vector of the track current position |
515 | Int_t copy ; |
7587f5a5 |
516 | |
bea63bea |
517 | Int_t tracknumber = gAlice->CurrentTrack() ; |
fa412d9b |
518 | Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() ); |
519 | TString name = fGeom->GetName() ; |
7eb9d12d |
520 | Int_t trackpid = gMC->TrackPid() ; |
fa412d9b |
521 | if ( name == "GPS2" ) { // ======> CPV is a GPS' PPSD |
522 | |
b2a60966 |
523 | if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell |
7587f5a5 |
524 | { |
525 | gMC->TrackPosition(pos) ; |
526 | xyze[0] = pos[0] ; |
527 | xyze[1] = pos[1] ; |
528 | xyze[2] = pos[2] ; |
bea63bea |
529 | xyze[3] = gMC->Edep() ; |
7587f5a5 |
530 | |
531 | if ( xyze[3] != 0 ) { // there is deposited energy |
532 | gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number |
533 | gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number |
534 | // 1-> Geom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper |
535 | // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower |
536 | gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell |
537 | gMC->CurrentVolID(relid[3]) ; // get the column number |
538 | |
539 | // get the absolute Id number |
540 | |
bea63bea |
541 | fGeom->RelToAbsNumbering(relid, absid) ; |
7587f5a5 |
542 | |
bea63bea |
543 | // add current hit to the hit list |
7eb9d12d |
544 | AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid); |
7587f5a5 |
545 | |
546 | } // there is deposited energy |
fa412d9b |
547 | } // We are inside the gas of the CPV |
548 | } // GPS2 configuration |
549 | |
550 | else if ( name == "IHEP" ) { // ======> CPV is a IHEP's one |
551 | |
552 | // Yuri Kharlov, 28 September 2000 |
553 | |
554 | if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") && |
555 | gMC->IsTrackEntering() && |
556 | gMC->TrackCharge() != 0) { |
557 | |
558 | // Charged track has just entered to the CPV sensitive plane |
559 | |
a3dfe79c |
560 | AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS"); |
fa412d9b |
561 | |
a3dfe79c |
562 | Int_t moduleNumber; |
563 | gMC->CurrentVolOffID(3,moduleNumber); |
564 | moduleNumber--; |
fa412d9b |
565 | |
566 | // Current position of the hit in the CPV module ref. system |
567 | |
568 | gMC -> TrackPosition(pos); |
569 | Float_t xyzm[3], xyzd[3], xyd[2]; |
cd461ab8 |
570 | Int_t i; |
571 | for (i=0; i<3; i++) xyzm[i] = pos[i]; |
fa412d9b |
572 | gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system |
573 | xyd[0] = xyzd[0]; |
574 | xyd[1] =-xyzd[2]; |
575 | |
576 | // Current momentum of the hit's track in the CPV module ref. system |
577 | |
578 | TLorentzVector pmom; |
579 | gMC -> TrackMomentum(pmom); |
580 | Float_t pm[3], pd[3]; |
cd461ab8 |
581 | for (i=0; i<3; i++) pm[i] = pmom[i]; |
fa412d9b |
582 | gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system |
583 | pmom[0] = pd[0]; |
584 | pmom[1] =-pd[1]; |
585 | pmom[2] =-pd[2]; |
586 | |
587 | // Current particle type of the hit's track |
588 | |
589 | Int_t ipart = gMC->TrackPid(); |
590 | |
591 | // Add the current particle in the list of the CPV hits. |
592 | |
d1b50469 |
593 | phos.GetCPVModule(moduleNumber).AddHit(fIshunt,tracknumber,pmom,xyd,ipart); |
fa412d9b |
594 | |
a3dfe79c |
595 | if (fDebugLevel == 1) { |
596 | printf("CPV hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n", |
597 | moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E()); |
598 | printf( " xy = (%8.4f, %8.4f) cm, ipart = %d\n", |
599 | xyd[0],xyd[1],ipart); |
600 | } |
fa412d9b |
601 | |
602 | // Digitize the current CPV hit: |
603 | |
604 | // 1. find pad response and |
605 | |
3d402178 |
606 | TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit |
a3dfe79c |
607 | CPVDigitize(pmom,xyd,moduleNumber,cpvDigits); |
fa412d9b |
608 | |
609 | Float_t xmean = 0; |
610 | Float_t zmean = 0; |
611 | Float_t qsum = 0; |
cd461ab8 |
612 | Int_t idigit,ndigits; |
fa412d9b |
613 | |
614 | // 2. go through the current digit list and sum digits in pads |
615 | |
616 | ndigits = cpvDigits->GetEntriesFast(); |
cd461ab8 |
617 | for (idigit=0; idigit<ndigits-1; idigit++) { |
3d402178 |
618 | AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit); |
fa412d9b |
619 | Float_t x1 = cpvDigit1->GetXpad() ; |
620 | Float_t z1 = cpvDigit1->GetYpad() ; |
621 | for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) { |
3d402178 |
622 | AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit); |
fa412d9b |
623 | Float_t x2 = cpvDigit2->GetXpad() ; |
624 | Float_t z2 = cpvDigit2->GetYpad() ; |
625 | if (x1==x2 && z1==z2) { |
626 | Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ; |
627 | cpvDigit2->SetQpad(qsum) ; |
628 | cpvDigits->RemoveAt(idigit) ; |
629 | } |
630 | } |
631 | } |
632 | cpvDigits->Compress() ; |
633 | |
634 | // 3. add digits to temporary hit list fTmpHits |
635 | |
636 | ndigits = cpvDigits->GetEntriesFast(); |
cd461ab8 |
637 | for (idigit=0; idigit<ndigits; idigit++) { |
3d402178 |
638 | AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit); |
a3dfe79c |
639 | relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number |
fa412d9b |
640 | relid[1] =-1 ; // means CPV |
641 | relid[2] = cpvDigit->GetXpad() ; // column number of a pad |
642 | relid[3] = cpvDigit->GetYpad() ; // row number of a pad |
643 | |
644 | // get the absolute Id number |
645 | fGeom->RelToAbsNumbering(relid, absid) ; |
646 | |
647 | // add current digit to the temporary hit list |
648 | xyze[0] = 0. ; |
649 | xyze[1] = 0. ; |
650 | xyze[2] = 0. ; |
651 | xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad |
652 | primary = -1; // No need in primary for CPV |
7eb9d12d |
653 | AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid); |
fa412d9b |
654 | |
655 | if (cpvDigit->GetQpad() > 0.02) { |
656 | xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5); |
657 | zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5); |
658 | qsum += cpvDigit->GetQpad(); |
659 | } |
660 | } |
661 | delete cpvDigits; |
662 | } |
663 | } // end of IHEP configuration |
7587f5a5 |
664 | |
fa412d9b |
665 | if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal |
666 | gMC->TrackPosition(pos) ; |
667 | xyze[0] = pos[0] ; |
668 | xyze[1] = pos[1] ; |
669 | xyze[2] = pos[2] ; |
670 | xyze[3] = gMC->Edep() ; |
671 | |
672 | if ( xyze[3] != 0 ) { |
673 | gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ; |
674 | relid[1] = 0 ; // means PBW04 |
675 | gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module |
676 | gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module |
677 | |
678 | // get the absolute Id number |
679 | |
680 | fGeom->RelToAbsNumbering(relid, absid) ; |
681 | |
682 | // add current hit to the hit list |
683 | |
7eb9d12d |
684 | AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid); |
fa412d9b |
685 | |
686 | } // there is deposited energy |
687 | } // we are inside a PHOS Xtal |
688 | } |
689 | |
690 | //____________________________________________________________________________ |
691 | void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits) |
692 | { |
693 | // ------------------------------------------------------------------------ |
694 | // Digitize one CPV hit: |
695 | // On input take exact 4-momentum p and position zxhit of the hit, |
696 | // find the pad response around this hit and |
697 | // put the amplitudes in the pads into array digits |
698 | // |
699 | // Author: Yuri Kharlov (after Serguei Sadovsky) |
700 | // 2 October 2000 |
701 | // ------------------------------------------------------------------------ |
702 | |
a3dfe79c |
703 | const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad) |
704 | const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e- |
705 | const Float_t kdEdx = 4.0; // Average energy loss in CPV; |
706 | const Int_t kNgamz = 5; // Ionization size in Z |
707 | const Int_t kNgamx = 9; // Ionization size in Phi |
708 | const Float_t kNoise = 0.03; // charge noise in one pad |
fa412d9b |
709 | |
710 | Float_t rnor1,rnor2; |
711 | |
712 | // Just a reminder on axes notation in the CPV module: |
713 | // axis Z goes along the beam |
714 | // axis X goes across the beam in the module plane |
715 | // axis Y is a normal to the module plane showing from the IP |
716 | |
717 | Float_t hitX = zxhit[0]; |
718 | Float_t hitZ =-zxhit[1]; |
719 | Float_t pX = p.Px(); |
720 | Float_t pZ =-p.Pz(); |
721 | Float_t pNorm = p.Py(); |
a3dfe79c |
722 | Float_t eloss = kdEdx; |
3d402178 |
723 | |
7eb9d12d |
724 | // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl; |
725 | |
fa412d9b |
726 | Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness(); |
727 | Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness(); |
728 | gRandom->Rannor(rnor1,rnor2); |
a3dfe79c |
729 | eloss *= (1 + kDetR*rnor1) * |
730 | TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2))); |
fa412d9b |
731 | Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2; |
732 | Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2; |
733 | Float_t zhit2 = zhit1 + dZY; |
734 | Float_t xhit2 = xhit1 + dXY; |
735 | |
a3dfe79c |
736 | Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in" |
737 | Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out" |
fa412d9b |
738 | |
739 | Int_t nIter; |
740 | Float_t zxe[3][5]; |
741 | if (iwht1==iwht2) { // incline 1-wire hit |
742 | nIter = 2; |
743 | zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2; |
a3dfe79c |
744 | zxe[1][0] = (iwht1 + 0.5) * kCelWr; |
745 | zxe[2][0] = eloss/2; |
fa412d9b |
746 | zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2; |
a3dfe79c |
747 | zxe[1][1] = (iwht1 + 0.5) * kCelWr; |
748 | zxe[2][1] = eloss/2; |
fa412d9b |
749 | } |
750 | else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit |
751 | nIter = 3; |
752 | Int_t iwht3 = (iwht1 + iwht2) / 2; |
a3dfe79c |
753 | Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1 |
754 | Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2 |
755 | Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3 |
fa412d9b |
756 | Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13 |
757 | Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23 |
758 | Float_t dxw1 = xhit1 - xwr13; |
759 | Float_t dxw2 = xhit2 - xwr23; |
a3dfe79c |
760 | Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr ); |
761 | Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr ); |
762 | Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr ); |
fa412d9b |
763 | zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2; |
764 | zxe[1][0] = xwht1; |
a3dfe79c |
765 | zxe[2][0] = eloss * egm1; |
fa412d9b |
766 | zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2; |
767 | zxe[1][1] = xwht2; |
a3dfe79c |
768 | zxe[2][1] = eloss * egm2; |
fa412d9b |
769 | zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1; |
770 | zxe[1][2] = xwht3; |
a3dfe79c |
771 | zxe[2][2] = eloss * egm3; |
fa412d9b |
772 | } |
773 | else { // incline 2-wire hit |
774 | nIter = 2; |
a3dfe79c |
775 | Float_t xwht1 = (iwht1 + 0.5) * kCelWr; |
776 | Float_t xwht2 = (iwht2 + 0.5) * kCelWr; |
fa412d9b |
777 | Float_t xwr12 = (xwht1 + xwht2) / 2; |
778 | Float_t dxw1 = xhit1 - xwr12; |
779 | Float_t dxw2 = xhit2 - xwr12; |
780 | Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) ); |
781 | Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) ); |
782 | zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2; |
783 | zxe[1][0] = xwht1; |
a3dfe79c |
784 | zxe[2][0] = eloss * egm1; |
fa412d9b |
785 | zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2; |
786 | zxe[1][1] = xwht2; |
a3dfe79c |
787 | zxe[2][1] = eloss * egm2; |
fa412d9b |
788 | } |
bea63bea |
789 | |
fa412d9b |
790 | // Finite size of ionization region |
791 | |
792 | Int_t nCellZ = fGeom->GetNumberOfPadsZ(); |
793 | Int_t nCellX = fGeom->GetNumberOfPadsPhi(); |
a3dfe79c |
794 | Int_t nz3 = (kNgamz+1)/2; |
795 | Int_t nx3 = (kNgamx+1)/2; |
796 | cpvDigits->Expand(nIter*kNgamx*kNgamz); |
fa412d9b |
797 | TClonesArray &ldigits = *(TClonesArray *)cpvDigits; |
798 | |
799 | for (Int_t iter=0; iter<nIter; iter++) { |
800 | |
801 | Float_t zhit = zxe[0][iter]; |
802 | Float_t xhit = zxe[1][iter]; |
803 | Float_t qhit = zxe[2][iter]; |
804 | Float_t zcell = zhit / fGeom->GetPadSizeZ(); |
805 | Float_t xcell = xhit / fGeom->GetPadSizePhi(); |
806 | if ( zcell<=0 || xcell<=0 || |
807 | zcell>=nCellZ || xcell>=nCellX) return; |
808 | Int_t izcell = (Int_t) zcell; |
809 | Int_t ixcell = (Int_t) xcell; |
810 | Float_t zc = zcell - izcell - 0.5; |
811 | Float_t xc = xcell - ixcell - 0.5; |
a3dfe79c |
812 | for (Int_t iz=1; iz<=kNgamz; iz++) { |
fa412d9b |
813 | Int_t kzg = izcell + iz - nz3; |
814 | if (kzg<=0 || kzg>nCellZ) continue; |
815 | Float_t zg = (Float_t)(iz-nz3) - zc; |
a3dfe79c |
816 | for (Int_t ix=1; ix<=kNgamx; ix++) { |
fa412d9b |
817 | Int_t kxg = ixcell + ix - nx3; |
818 | if (kxg<=0 || kxg>nCellX) continue; |
819 | Float_t xg = (Float_t)(ix-nx3) - xc; |
820 | |
821 | // Now calculate pad response |
822 | Float_t qpad = CPVPadResponseFunction(qhit,zg,xg); |
a3dfe79c |
823 | qpad += kNoise*rnor2; |
fa412d9b |
824 | if (qpad<0) continue; |
825 | |
826 | // Fill the array with pad response ID and amplitude |
3d402178 |
827 | new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad); |
fa412d9b |
828 | } |
fa412d9b |
829 | } |
fa412d9b |
830 | } |
831 | } |
832 | |
833 | //____________________________________________________________________________ |
834 | Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) { |
835 | // ------------------------------------------------------------------------ |
836 | // Calculate the amplitude in one CPV pad using the |
837 | // cumulative pad response function |
838 | // Author: Yuri Kharlov (after Serguei Sadovski) |
839 | // 3 October 2000 |
840 | // ------------------------------------------------------------------------ |
841 | |
842 | Double_t dz = fGeom->GetPadSizeZ() / 2; |
843 | Double_t dx = fGeom->GetPadSizePhi() / 2; |
844 | Double_t z = zhit * fGeom->GetPadSizeZ(); |
845 | Double_t x = xhit * fGeom->GetPadSizePhi(); |
846 | Double_t amplitude = qhit * |
847 | (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) - |
848 | CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx)); |
849 | return (Float_t)amplitude; |
7587f5a5 |
850 | } |
851 | |
fa412d9b |
852 | //____________________________________________________________________________ |
853 | Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) { |
854 | // ------------------------------------------------------------------------ |
855 | // Cumulative pad response function |
856 | // It includes several terms from the CF decomposition in electrostatics |
857 | // Note: this cumulative function is wrong since omits some terms |
858 | // but the cell amplitude obtained with it is correct because |
859 | // these omitting terms cancel |
860 | // Author: Yuri Kharlov (after Serguei Sadovski) |
861 | // 3 October 2000 |
862 | // ------------------------------------------------------------------------ |
863 | |
a3dfe79c |
864 | const Double_t kA=1.0; |
865 | const Double_t kB=0.7; |
fa412d9b |
866 | |
867 | Double_t r2 = x*x + y*y; |
868 | Double_t xy = x*y; |
869 | Double_t cumulPRF = 0; |
870 | for (Int_t i=0; i<=4; i++) { |
a3dfe79c |
871 | Double_t b1 = (2*i + 1) * kB; |
fa412d9b |
872 | cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) ); |
873 | } |
a3dfe79c |
874 | cumulPRF *= kA/(2*TMath::Pi()); |
fa412d9b |
875 | return cumulPRF; |
876 | } |
7eb9d12d |
877 | |