]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSv1.cxx
Added virtual IsVersion. Print of version defered to Init in AliPHOS
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv1.cxx
CommitLineData
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
20// Layout EMC + PPSD has name GPS2
21// Produces cumulated hits (no hits) and digits
22//*-- Author: Yves Schutz (SUBATECH)
b2a60966 23
7587f5a5 24
25// --- ROOT system ---
bea63bea 26
27#include "TBRIK.h"
28#include "TNode.h"
7587f5a5 29#include "TRandom.h"
30
5f20d3fb 31
7587f5a5 32// --- Standard library ---
33
de9ec31b 34#include <stdio.h>
35#include <string.h>
36#include <stdlib.h>
37#include <strstream.h>
7587f5a5 38
39// --- AliRoot header files ---
40
41#include "AliPHOSv1.h"
42#include "AliPHOSHit.h"
43#include "AliPHOSDigit.h"
bea63bea 44#include "AliPHOSReconstructioner.h"
7587f5a5 45#include "AliRun.h"
46#include "AliConst.h"
47
48ClassImp(AliPHOSv1)
49
bea63bea 50//____________________________________________________________________________
51AliPHOSv1::AliPHOSv1()
52{
5f20d3fb 53 // ctor
bea63bea 54 fNTmpHits = 0 ;
55 fTmpHits = 0 ;
56}
57
7587f5a5 58//____________________________________________________________________________
59AliPHOSv1::AliPHOSv1(const char *name, const char *title):
e04976bd 60AliPHOSv0(name,title)
7587f5a5 61{
5f20d3fb 62 // ctor : title is used to identify the layout
63 // GPS2 = 5 modules (EMC + PPSD)
64 // We use 2 arrays of hits :
65 //
66 // - fHits (the "normal" one), which retains the hits associated with
67 // the current primary particle being tracked
68 // (this array is reset after each primary has been tracked).
69 //
70 // - fTmpHits, which retains all the hits of the current event. It
71 // is used for the digitization part.
4a2ca5e9 72
5f20d3fb 73 fPinElectronicNoise = 0.010 ;
b95514c0 74 fDigitThreshold = 0.1 ; // 1 GeV
5f20d3fb 75
76 // We do not want to save in TreeH the raw hits
77 // But save the cumulated hits instead (need to create the branch myself)
78 // It is put in the Digit Tree because the TreeH is filled after each primary
79 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
80
81 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
82
83 fNTmpHits = fNhits = 0 ;
84
85 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
86
87
88 fIshunt = 1 ; // All hits are associated with primary particles
e04976bd 89
5f20d3fb 90}
91
92//____________________________________________________________________________
93AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
94 AliPHOSv0(name,title)
95{
96 // ctor : title is used to identify the layout
97 // GPS2 = 5 modules (EMC + PPSD)
98 // We use 2 arrays of hits :
99 //
100 // - fHits (the "normal" one), which retains the hits associated with
101 // the current primary particle being tracked
102 // (this array is reset after each primary has been tracked).
103 //
104 // - fTmpHits, which retains all the hits of the current event. It
105 // is used for the digitization part.
106
107 fPinElectronicNoise = 0.010 ;
108
109 // We do not want to save in TreeH the raw hits
110 //fHits = new TClonesArray("AliPHOSHit",100) ;
111
112 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
113 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
114
115 fNTmpHits = fNhits = 0 ;
116
117 fIshunt = 1 ; // All hits are associated with primary particles
118
119 // gets an instance of the geometry parameters class
120 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
121
122 if (fGeom->IsInitialized() )
88bdfa12 123 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
5f20d3fb 124 else
88bdfa12 125 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
5f20d3fb 126
127 // Defining the PHOS Reconstructioner
128
129 fReconstructioner = Reconstructioner ;
130
7587f5a5 131}
b2a60966 132
7587f5a5 133//____________________________________________________________________________
bea63bea 134AliPHOSv1::~AliPHOSv1()
b2a60966 135{
bea63bea 136 // dtor
5f20d3fb 137
8dfa469d 138 if ( fTmpHits) {
139 fTmpHits->Delete() ;
140 delete fTmpHits ;
141 fTmpHits = 0 ;
142 }
5f20d3fb 143
144 if ( fEmcRecPoints ) {
8dfa469d 145 fEmcRecPoints->Delete() ;
146 delete fEmcRecPoints ;
147 fEmcRecPoints = 0 ;
148 }
5f20d3fb 149
8dfa469d 150 if ( fPpsdRecPoints ) {
151 fPpsdRecPoints->Delete() ;
152 delete fPpsdRecPoints ;
153 fPpsdRecPoints = 0 ;
154 }
5f20d3fb 155
8dfa469d 156 if ( fTrackSegments ) {
157 fTrackSegments->Delete() ;
158 delete fTrackSegments ;
159 fTrackSegments = 0 ;
160 }
5f20d3fb 161
7587f5a5 162}
163
7587f5a5 164//____________________________________________________________________________
bea63bea 165void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
166{
167 // Add a hit to the hit list.
5f20d3fb 168 // A PHOS hit is the sum of all hits in a single crystal
169 // or in a single PPSD gas cell
bea63bea 170
5f20d3fb 171 Int_t hitCounter ;
172 TClonesArray &ltmphits = *fTmpHits ;
bea63bea 173 AliPHOSHit *newHit ;
5f20d3fb 174 AliPHOSHit *curHit ;
175 Bool_t deja = kFALSE ;
bea63bea 176
5f20d3fb 177 // In any case, fills the fTmpHit TClonesArray (with "accumulated hits")
bea63bea 178
179 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
180
5f20d3fb 181 // We do not want to save in TreeH the raw hits
bea63bea 182 // TClonesArray &lhits = *fHits;
bea63bea 183
5f20d3fb 184 for ( hitCounter = 0 ; hitCounter < fNTmpHits && !deja ; hitCounter++ ) {
185 curHit = (AliPHOSHit*) ltmphits[hitCounter] ;
186 if( *curHit == *newHit ) {
187 *curHit = *curHit + *newHit ;
188 deja = kTRUE ;
189 }
190 }
191
192 if ( !deja ) {
193 new(ltmphits[fNTmpHits]) AliPHOSHit(*newHit) ;
194 fNTmpHits++ ;
195 }
196
bea63bea 197 // We do not want to save in TreeH the raw hits
198 // new(lhits[fNhits]) AliPHOSHit(*newHit) ;
199 // fNhits++ ;
200
201 // Please note that the fTmpHits array must survive up to the
202 // end of the events, so it does not appear e.g. in ResetHits() (
203 // which is called at the end of each primary).
204
205 delete newHit;
206
207}
208
5f20d3fb 209//___________________________________________________________________________
210Int_t AliPHOSv1::Digitize(Float_t Energy)
211{
212 // Applies the energy calibration
213
214 Float_t fB = 100000000. ;
215 Float_t fA = 0. ;
216 Int_t chan = Int_t(fA + Energy*fB ) ;
217 return chan ;
218}
bea63bea 219
220//___________________________________________________________________________
221void AliPHOSv1::FinishEvent()
222{
223 // Makes the digits from the sum of summed hit in a single crystal or PPSD gas cell
224 // Adds to the energy the electronic noise
225 // Keeps digits with energy above fDigitThreshold
226
227 // Save the cumulated hits instead of raw hits (need to create the branch myself)
228 // It is put in the Digit Tree because the TreeH is filled after each primary
229 // and the TreeD at the end of the event.
230
5f20d3fb 231
bea63bea 232 Int_t i ;
233 Int_t relid[4];
234 Int_t j ;
235 TClonesArray &lDigits = *fDigits ;
236 AliPHOSHit * hit ;
237 AliPHOSDigit * newdigit ;
238 AliPHOSDigit * curdigit ;
239 Bool_t deja = kFALSE ;
240
5f20d3fb 241 for ( i = 0 ; i < fNTmpHits ; i++ ) {
242 hit = (AliPHOSHit*)fTmpHits->At(i) ;
4a2ca5e9 243
244 // Assign primary number only if contribution is significant
245 if( hit->GetEnergy() > fDigitThreshold)
246 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
247 else
248 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
bea63bea 249 deja =kFALSE ;
250 for ( j = 0 ; j < fNdigits ; j++) {
251 curdigit = (AliPHOSDigit*) lDigits[j] ;
252 if ( *curdigit == *newdigit) {
253 *curdigit = *curdigit + *newdigit ;
254 deja = kTRUE ;
255 }
256 }
257 if ( !deja ) {
258 new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ;
259 fNdigits++ ;
260 }
261
262 delete newdigit ;
263 }
264
265 // Noise induced by the PIN diode of the PbWO crystals
266
267 Float_t energyandnoise ;
268 for ( i = 0 ; i < fNdigits ; i++ ) {
269 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
270 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
271
272 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
273 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
274
275 if (energyandnoise < 0 )
276 energyandnoise = 0 ;
277
278 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
279 fDigits->RemoveAt(i) ;
280 }
281 }
282
283 fDigits->Compress() ;
284
285 fNdigits = fDigits->GetEntries() ;
286 for (i = 0 ; i < fNdigits ; i++) {
287 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
288 newdigit->SetIndexInList(i) ;
289 }
432ec55f 290
5f20d3fb 291}
292
293//___________________________________________________________________________
294void AliPHOSv1::MakeBranch(Option_t* opt)
295{
296 // Create new branche in the current Root Tree in the digit Tree
5f20d3fb 297 AliDetector::MakeBranch(opt) ;
298
299 char branchname[10];
300 sprintf(branchname,"%s",GetName());
301 char *cdD = strstr(opt,"D");
302 if (fDigits && gAlice->TreeD() && cdD) {
303 gAlice->TreeD()->Branch(branchname, &fDigits, fBufferSize);
304 }
305
306 // Create new branche PHOSCH in the current Root Tree in the digit Tree for accumulated Hits
307 if ( ! (gAlice->IsLegoRun()) ) { // only when not in lego plot mode
308 if ( fTmpHits && gAlice->TreeD() && cdD) {
309 char branchname[10] ;
310 sprintf(branchname, "%sCH", GetName()) ;
311 gAlice->TreeD()->Branch(branchname, &fTmpHits, fBufferSize) ;
312 }
313 }
bea63bea 314
315}
316
5f20d3fb 317//_____________________________________________________________________________
318void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
319{
320 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
321 // 2. Creates TreeR with a branch for each list
322 // 3. Steers the reconstruction processes
323 // 4. Saves the 3 lists in TreeR
324 // 5. Write the Tree to File
325
326 fReconstructioner = Reconstructioner ;
327
328 char branchname[10] ;
329
330 // 1.
331
332 // gAlice->MakeTree("R") ;
333 Int_t splitlevel = 0 ;
334
335 if (fEmcRecPoints) {
336 fEmcRecPoints->Delete() ;
337 delete fEmcRecPoints ;
338 fEmcRecPoints = 0 ;
339 }
340
341 // fEmcRecPoints= new AliPHOSRecPoint::RecPointsList("AliPHOSEmcRecPoint", 1000) ; if TClonesArray
342 fEmcRecPoints= new AliPHOSRecPoint::RecPointsList(2000) ;
343
344 if ( fEmcRecPoints && gAlice->TreeR() ) {
345 sprintf(branchname,"%sEmcRP",GetName()) ;
346
347 // gAlice->TreeR()->Branch(branchname, &fEmcRecPoints, fBufferSize); if TClonesArray
348 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
349 }
350
351 if (fPpsdRecPoints) {
352 fPpsdRecPoints->Delete() ;
353 delete fPpsdRecPoints ;
354 fPpsdRecPoints = 0 ;
355 }
356
357 // fPpsdRecPoints = new AliPHOSRecPoint::RecPointsList("AliPHOSPpsdRecPoint", 1000) ; if TClonesArray
358 fPpsdRecPoints = new AliPHOSRecPoint::RecPointsList(2000) ;
359
360 if ( fPpsdRecPoints && gAlice->TreeR() ) {
361 sprintf(branchname,"%sPpsdRP",GetName()) ;
362
363 // gAlice->TreeR()->Branch(branchname, &fPpsdRecPoints, fBufferSize); if TClonesArray
364 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
365 }
366
367 if (fTrackSegments) {
368 fTrackSegments->Delete() ;
369 delete fTrackSegments ;
370 fTrackSegments = 0 ;
371 }
372
373 fTrackSegments = new AliPHOSTrackSegment::TrackSegmentsList("AliPHOSTrackSegment", 2000) ;
374 if ( fTrackSegments && gAlice->TreeR() ) {
375 sprintf(branchname,"%sTS",GetName()) ;
376 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
377 }
378
379 if (fRecParticles) {
380 fRecParticles->Delete() ;
381 delete fRecParticles ;
382 fRecParticles = 0 ;
383 }
384 fRecParticles = new AliPHOSRecParticle::RecParticlesList("AliPHOSRecParticle", 2000) ;
385 if ( fRecParticles && gAlice->TreeR() ) {
386 sprintf(branchname,"%sRP",GetName()) ;
387 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
388 }
389
390 // 3.
391
392 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
393
394 // 4. Expand or Shrink the arrays to the proper size
395
396 Int_t size ;
397
398 size = fEmcRecPoints->GetEntries() ;
399 fEmcRecPoints->Expand(size) ;
400
401 size = fPpsdRecPoints->GetEntries() ;
402 fPpsdRecPoints->Expand(size) ;
403
404 size = fTrackSegments->GetEntries() ;
405 fTrackSegments->Expand(size) ;
406
407 size = fRecParticles->GetEntries() ;
408 fRecParticles->Expand(size) ;
409
410 gAlice->TreeR()->Fill() ;
5f20d3fb 411 // 5.
412
4a2ca5e9 413 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
5f20d3fb 414
415 // Deleting reconstructed objects
416 ResetReconstruction();
417
418}
419
420//____________________________________________________________________________
421void AliPHOSv1::ResetDigits()
422{
423 // May sound strange, but cumulative hits are store in digits Tree
424 AliDetector::ResetDigits();
425 if( fTmpHits ) {
426 fTmpHits->Delete();
427 fNTmpHits = 0 ;
428 }
429}
430//____________________________________________________________________________
431void AliPHOSv1::ResetReconstruction()
432{
433 // Deleting reconstructed objects
434
435 if ( fEmcRecPoints ) fEmcRecPoints->Delete();
436 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
437 if ( fTrackSegments ) fTrackSegments->Delete();
438 if ( fRecParticles ) fRecParticles->Delete();
439
440}
5f20d3fb 441
442//____________________________________________________________________________
443void AliPHOSv1::SetTreeAddress()
444{
445 // TBranch *branch;
446 AliPHOS::SetTreeAddress();
447
448 // //Branch address for TreeR: RecPpsdRecPoint
449// TTree *treeR = gAlice->TreeR();
450// if ( treeR && fPpsdRecPoints ) {
451// branch = treeR->GetBranch("PHOSPpsdRP");
452// if (branch) branch->SetAddress(&fPpsdRecPoints) ;
453 // }
454}
455
456//____________________________________________________________________________
457
7587f5a5 458void AliPHOSv1::StepManager(void)
459{
b2a60966 460 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
b2a60966 461
7587f5a5 462 Int_t relid[4] ; // (box, layer, row, column) indices
463 Float_t xyze[4] ; // position wrt MRS and energy deposited
464 TLorentzVector pos ;
bea63bea 465 Int_t copy ;
7587f5a5 466
bea63bea 467 Int_t tracknumber = gAlice->CurrentTrack() ;
5f20d3fb 468 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
bea63bea 469 TString name = fGeom->GetName() ;
7587f5a5 470 if ( name == "GPS2" ) { // the CPV is a PPSD
b2a60966 471 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
7587f5a5 472 {
473 gMC->TrackPosition(pos) ;
474 xyze[0] = pos[0] ;
475 xyze[1] = pos[1] ;
476 xyze[2] = pos[2] ;
bea63bea 477 xyze[3] = gMC->Edep() ;
7587f5a5 478
479 if ( xyze[3] != 0 ) { // there is deposited energy
480 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
481 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
482 // 1-> Geom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
483 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
484 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
485 gMC->CurrentVolID(relid[3]) ; // get the column number
486
487 // get the absolute Id number
488
489 Int_t absid ;
bea63bea 490 fGeom->RelToAbsNumbering(relid, absid) ;
7587f5a5 491
bea63bea 492 // add current hit to the hit list
493 AddHit(fIshunt, primary, tracknumber, absid, xyze);
7587f5a5 494
495 } // there is deposited energy
496 } // We are inside the gas of the CPV
497 } // GPS2 configuration
498
bea63bea 499 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) // We are inside a PBWO crystal
7587f5a5 500 {
501 gMC->TrackPosition(pos) ;
502 xyze[0] = pos[0] ;
503 xyze[1] = pos[1] ;
504 xyze[2] = pos[2] ;
a333c916 505 xyze[3] = gMC->Edep() ;
7587f5a5 506
7587f5a5 507 if ( xyze[3] != 0 ) {
508 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
bea63bea 509 relid[1] = 0 ; // means PBW04
7587f5a5 510 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
511 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
512
5f20d3fb 513 // get the absolute Id number
514
7587f5a5 515 Int_t absid ;
bea63bea 516 fGeom->RelToAbsNumbering(relid, absid) ;
517
5f20d3fb 518 // add current hit to the hit list
bea63bea 519
5f20d3fb 520 AddHit(fIshunt, primary,tracknumber, absid, xyze);
521
7587f5a5 522 } // there is deposited energy
5f20d3fb 523 } // we are inside a PHOS Xtal
7587f5a5 524}
525