1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
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)
25 // --- ROOT system ---
33 // --- Standard library ---
38 #include <strstream.h>
40 // --- AliRoot header files ---
42 #include "AliPHOSv1.h"
43 #include "AliPHOSHit.h"
44 #include "AliPHOSDigit.h"
45 #include "AliPHOSReconstructioner.h"
52 //____________________________________________________________________________
53 AliPHOSv1::AliPHOSv1()
62 //____________________________________________________________________________
63 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
66 // ctor : title is used to identify the layout
67 // GPS2 = 5 modules (EMC + PPSD)
68 // We use 2 arrays of hits :
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).
74 // - fTmpHits, which retains all the hits of the current event. It
75 // is used for the digitization part.
77 fPinElectronicNoise = 0.010 ;
78 fDigitThreshold = 0.1 ; // 1 GeV
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() ).
85 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
87 fNTmpHits = fNhits = 0 ;
89 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
92 fIshunt = 1 ; // All hits are associated with primary particles
96 //____________________________________________________________________________
97 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
100 // ctor : title is used to identify the layout
101 // GPS2 = 5 modules (EMC + PPSD)
102 // We use 2 arrays of hits :
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).
108 // - fTmpHits, which retains all the hits of the current event. It
109 // is used for the digitization part.
111 fPinElectronicNoise = 0.010 ;
113 // We do not want to save in TreeH the raw hits
114 //fHits = new TClonesArray("AliPHOSHit",100) ;
116 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
117 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
119 fNTmpHits = fNhits = 0 ;
121 fIshunt = 1 ; // All hits are associated with primary particles
123 // gets an instance of the geometry parameters class
124 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
126 if (fGeom->IsInitialized() )
127 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
129 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
131 // Defining the PHOS Reconstructioner
133 fReconstructioner = Reconstructioner ;
137 //____________________________________________________________________________
138 AliPHOSv1::~AliPHOSv1()
148 // if ( fEmcRecPoints ) {
149 // fEmcRecPoints->Delete() ;
150 // delete fEmcRecPoints ;
151 // fEmcRecPoints = 0 ;
154 // if ( fPpsdRecPoints ) {
155 // fPpsdRecPoints->Delete() ;
156 // delete fPpsdRecPoints ;
157 // fPpsdRecPoints = 0 ;
160 // if ( fTrackSegments ) {
161 // fTrackSegments->Delete() ;
162 // delete fTrackSegments ;
163 // fTrackSegments = 0 ;
168 //____________________________________________________________________________
169 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
171 // Add a hit to the hit list.
172 // A PHOS hit is the sum of all hits in a single crystal
173 // or in a single PPSD gas cell
176 TClonesArray <mphits = *fTmpHits ;
179 Bool_t deja = kFALSE ;
181 // In any case, fills the fTmpHit TClonesArray (with "accumulated hits")
183 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
185 // We do not want to save in TreeH the raw hits
186 // TClonesArray &lhits = *fHits;
188 for ( hitCounter = 0 ; hitCounter < fNTmpHits && !deja ; hitCounter++ ) {
189 curHit = (AliPHOSHit*) ltmphits[hitCounter] ;
190 if( *curHit == *newHit ) {
191 *curHit = *curHit + *newHit ;
197 new(ltmphits[fNTmpHits]) AliPHOSHit(*newHit) ;
201 // We do not want to save in TreeH the raw hits
202 // new(lhits[fNhits]) AliPHOSHit(*newHit) ;
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).
213 //___________________________________________________________________________
214 Int_t AliPHOSv1::Digitize(Float_t Energy)
216 // Applies the energy calibration
218 Float_t fB = 100000000. ;
220 Int_t chan = Int_t(fA + Energy*fB ) ;
224 //___________________________________________________________________________
225 void AliPHOSv1::FinishEvent()
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
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.
239 TClonesArray &lDigits = *fDigits ;
241 AliPHOSDigit * newdigit ;
242 AliPHOSDigit * curdigit ;
243 Bool_t deja = kFALSE ;
245 for ( i = 0 ; i < fNTmpHits ; i++ ) {
246 hit = (AliPHOSHit*)fTmpHits->At(i) ;
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() ) ) ;
252 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
254 for ( j = 0 ; j < fNdigits ; j++) {
255 curdigit = (AliPHOSDigit*) lDigits[j] ;
256 if ( *curdigit == *newdigit) {
257 *curdigit = *curdigit + *newdigit ;
262 new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ;
269 // Noise induced by the PIN diode of the PbWO crystals
271 Float_t energyandnoise ;
272 for ( i = 0 ; i < fNdigits ; i++ ) {
273 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
274 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
276 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
277 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
279 if (energyandnoise < 0 )
282 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
283 fDigits->RemoveAt(i) ;
287 fDigits->Compress() ;
289 fNdigits = fDigits->GetEntries() ;
290 for (i = 0 ; i < fNdigits ; i++) {
291 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
292 newdigit->SetIndexInList(i) ;
297 //___________________________________________________________________________
298 void AliPHOSv1::MakeBranch(Option_t* opt)
300 // Create new branche in the current Root Tree in the digit Tree
301 AliDetector::MakeBranch(opt) ;
304 sprintf(branchname,"%s",GetName());
305 char *cdD = strstr(opt,"D");
306 if (fDigits && gAlice->TreeD() && cdD) {
307 gAlice->TreeD()->Branch(branchname, &fDigits, fBufferSize);
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) ;
321 //_____________________________________________________________________________
322 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
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
330 fReconstructioner = Reconstructioner ;
332 char branchname[10] ;
336 // gAlice->MakeTree("R") ;
337 Int_t splitlevel = 0 ;
339 fEmcRecPoints->Delete() ;
341 if ( fEmcRecPoints && gAlice->TreeR() ) {
342 sprintf(branchname,"%sEmcRP",GetName()) ;
344 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
347 fPpsdRecPoints->Delete() ;
350 if ( fPpsdRecPoints && gAlice->TreeR() ) {
351 sprintf(branchname,"%sPpsdRP",GetName()) ;
353 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
356 fTrackSegments->Delete() ;
358 if ( fTrackSegments && gAlice->TreeR() ) {
359 sprintf(branchname,"%sTS",GetName()) ;
360 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
363 fRecParticles->Delete() ;
365 if ( fRecParticles && gAlice->TreeR() ) {
366 sprintf(branchname,"%sRP",GetName()) ;
367 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
372 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
374 // 4. Expand or Shrink the arrays to the proper size
378 size = fEmcRecPoints->GetEntries() ;
379 fEmcRecPoints->Expand(size) ;
381 size = fPpsdRecPoints->GetEntries() ;
382 fPpsdRecPoints->Expand(size) ;
384 size = fTrackSegments->GetEntries() ;
385 fTrackSegments->Expand(size) ;
387 size = fRecParticles->GetEntries() ;
388 fRecParticles->Expand(size) ;
390 gAlice->TreeR()->Fill() ;
393 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
395 // Deleting reconstructed objects
396 ResetReconstruction();
400 //____________________________________________________________________________
401 void AliPHOSv1::ResetDigits()
403 // May sound strange, but cumulative hits are store in digits Tree
404 AliDetector::ResetDigits();
410 //____________________________________________________________________________
411 void AliPHOSv1::ResetReconstruction()
413 // Deleting reconstructed objects
415 if ( fEmcRecPoints ) fEmcRecPoints->Delete();
416 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
417 if ( fTrackSegments ) fTrackSegments->Delete();
418 if ( fRecParticles ) fRecParticles->Delete();
422 //____________________________________________________________________________
423 void AliPHOSv1::SetTreeAddress()
426 AliPHOS::SetTreeAddress();
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) ;
436 //____________________________________________________________________________
438 void AliPHOSv1::StepManager(void)
440 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
442 Int_t relid[4] ; // (box, layer, row, column) indices
443 Float_t xyze[4] ; // position wrt MRS and energy deposited
447 Int_t tracknumber = gAlice->CurrentTrack() ;
448 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
449 TString name = fGeom->GetName() ;
450 if ( name == "GPS2" ) { // the CPV is a PPSD
451 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
453 gMC->TrackPosition(pos) ;
457 xyze[3] = gMC->Edep() ;
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
467 // get the absolute Id number
470 fGeom->RelToAbsNumbering(relid, absid) ;
472 // add current hit to the hit list
473 AddHit(fIshunt, primary, tracknumber, absid, xyze);
475 } // there is deposited energy
476 } // We are inside the gas of the CPV
477 } // GPS2 configuration
479 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) // We are inside a PBWO crystal
481 gMC->TrackPosition(pos) ;
485 xyze[3] = gMC->Edep() ;
487 if ( xyze[3] != 0 ) {
488 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
489 relid[1] = 0 ; // means PBW04
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
493 // get the absolute Id number
496 fGeom->RelToAbsNumbering(relid, absid) ;
498 // add current hit to the hit list
500 AddHit(fIshunt, primary,tracknumber, absid, xyze);
502 } // there is deposited energy
503 } // we are inside a PHOS Xtal