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