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 ---
32 // --- Standard library ---
37 #include <strstream.h>
39 // --- AliRoot header files ---
41 #include "AliPHOSv1.h"
42 #include "AliPHOSHit.h"
43 #include "AliPHOSDigit.h"
44 #include "AliPHOSReconstructioner.h"
50 //____________________________________________________________________________
51 AliPHOSv1::AliPHOSv1()
60 //____________________________________________________________________________
61 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
64 // ctor : title is used to identify the layout
65 // GPS2 = 5 modules (EMC + PPSD)
66 // We use 2 arrays of hits :
68 // - fHits (the "normal" one), which retains the hits associated with
69 // the current primary particle being tracked
70 // (this array is reset after each primary has been tracked).
72 // - fTmpHits, which retains all the hits of the current event. It
73 // is used for the digitization part.
75 fPinElectronicNoise = 0.010 ;
76 fDigitThreshold = 0.1 ; // 1 GeV
78 // We do not want to save in TreeH the raw hits
79 // But save the cumulated hits instead (need to create the branch myself)
80 // It is put in the Digit Tree because the TreeH is filled after each primary
81 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
83 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
85 fNTmpHits = fNhits = 0 ;
87 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
90 fIshunt = 1 ; // All hits are associated with primary particles
94 //____________________________________________________________________________
95 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
98 // ctor : title is used to identify the layout
99 // GPS2 = 5 modules (EMC + PPSD)
100 // We use 2 arrays of hits :
102 // - fHits (the "normal" one), which retains the hits associated with
103 // the current primary particle being tracked
104 // (this array is reset after each primary has been tracked).
106 // - fTmpHits, which retains all the hits of the current event. It
107 // is used for the digitization part.
109 fPinElectronicNoise = 0.010 ;
111 // We do not want to save in TreeH the raw hits
112 //fHits = new TClonesArray("AliPHOSHit",100) ;
114 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
115 fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
117 fNTmpHits = fNhits = 0 ;
119 fIshunt = 1 ; // All hits are associated with primary particles
121 // gets an instance of the geometry parameters class
122 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
124 if (fGeom->IsInitialized() )
125 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
127 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
129 // Defining the PHOS Reconstructioner
131 fReconstructioner = Reconstructioner ;
135 //____________________________________________________________________________
136 AliPHOSv1::~AliPHOSv1()
146 // if ( fEmcRecPoints ) {
147 // fEmcRecPoints->Delete() ;
148 // delete fEmcRecPoints ;
149 // fEmcRecPoints = 0 ;
152 // if ( fPpsdRecPoints ) {
153 // fPpsdRecPoints->Delete() ;
154 // delete fPpsdRecPoints ;
155 // fPpsdRecPoints = 0 ;
158 // if ( fTrackSegments ) {
159 // fTrackSegments->Delete() ;
160 // delete fTrackSegments ;
161 // fTrackSegments = 0 ;
166 //____________________________________________________________________________
167 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
169 // Add a hit to the hit list.
170 // A PHOS hit is the sum of all hits in a single crystal
171 // or in a single PPSD gas cell
174 TClonesArray <mphits = *fTmpHits ;
177 Bool_t deja = kFALSE ;
179 // In any case, fills the fTmpHit TClonesArray (with "accumulated hits")
181 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
183 // We do not want to save in TreeH the raw hits
184 // TClonesArray &lhits = *fHits;
186 for ( hitCounter = 0 ; hitCounter < fNTmpHits && !deja ; hitCounter++ ) {
187 curHit = (AliPHOSHit*) ltmphits[hitCounter] ;
188 if( *curHit == *newHit ) {
189 *curHit = *curHit + *newHit ;
195 new(ltmphits[fNTmpHits]) AliPHOSHit(*newHit) ;
199 // We do not want to save in TreeH the raw hits
200 // new(lhits[fNhits]) AliPHOSHit(*newHit) ;
203 // Please note that the fTmpHits array must survive up to the
204 // end of the events, so it does not appear e.g. in ResetHits() (
205 // which is called at the end of each primary).
211 //___________________________________________________________________________
212 Int_t AliPHOSv1::Digitize(Float_t Energy)
214 // Applies the energy calibration
216 Float_t fB = 100000000. ;
218 Int_t chan = Int_t(fA + Energy*fB ) ;
222 //___________________________________________________________________________
223 void AliPHOSv1::FinishEvent()
225 // Makes the digits from the sum of summed hit in a single crystal or PPSD gas cell
226 // Adds to the energy the electronic noise
227 // Keeps digits with energy above fDigitThreshold
229 // Save the cumulated hits instead of raw hits (need to create the branch myself)
230 // It is put in the Digit Tree because the TreeH is filled after each primary
231 // and the TreeD at the end of the event.
237 TClonesArray &lDigits = *fDigits ;
239 AliPHOSDigit * newdigit ;
240 AliPHOSDigit * curdigit ;
241 Bool_t deja = kFALSE ;
243 for ( i = 0 ; i < fNTmpHits ; i++ ) {
244 hit = (AliPHOSHit*)fTmpHits->At(i) ;
246 // Assign primary number only if contribution is significant
247 if( hit->GetEnergy() > fDigitThreshold)
248 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
250 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
252 for ( j = 0 ; j < fNdigits ; j++) {
253 curdigit = (AliPHOSDigit*) lDigits[j] ;
254 if ( *curdigit == *newdigit) {
255 *curdigit = *curdigit + *newdigit ;
260 new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ;
267 // Noise induced by the PIN diode of the PbWO crystals
269 Float_t energyandnoise ;
270 for ( i = 0 ; i < fNdigits ; i++ ) {
271 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
272 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
274 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
275 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
277 if (energyandnoise < 0 )
280 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
281 fDigits->RemoveAt(i) ;
285 fDigits->Compress() ;
287 fNdigits = fDigits->GetEntries() ;
288 for (i = 0 ; i < fNdigits ; i++) {
289 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
290 newdigit->SetIndexInList(i) ;
295 //___________________________________________________________________________
296 void AliPHOSv1::MakeBranch(Option_t* opt)
298 // Create new branche in the current Root Tree in the digit Tree
299 AliDetector::MakeBranch(opt) ;
302 sprintf(branchname,"%s",GetName());
303 char *cdD = strstr(opt,"D");
304 if (fDigits && gAlice->TreeD() && cdD) {
305 gAlice->TreeD()->Branch(branchname, &fDigits, fBufferSize);
308 // Create new branche PHOSCH in the current Root Tree in the digit Tree for accumulated Hits
309 if ( ! (gAlice->IsLegoRun()) ) { // only when not in lego plot mode
310 if ( fTmpHits && gAlice->TreeD() && cdD) {
311 char branchname[10] ;
312 sprintf(branchname, "%sCH", GetName()) ;
313 gAlice->TreeD()->Branch(branchname, &fTmpHits, fBufferSize) ;
319 //_____________________________________________________________________________
320 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
322 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
323 // 2. Creates TreeR with a branch for each list
324 // 3. Steers the reconstruction processes
325 // 4. Saves the 3 lists in TreeR
326 // 5. Write the Tree to File
328 fReconstructioner = Reconstructioner ;
330 char branchname[10] ;
334 // gAlice->MakeTree("R") ;
335 Int_t splitlevel = 0 ;
337 fEmcRecPoints->Delete() ;
339 if ( fEmcRecPoints && gAlice->TreeR() ) {
340 sprintf(branchname,"%sEmcRP",GetName()) ;
342 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
345 fPpsdRecPoints->Delete() ;
348 if ( fPpsdRecPoints && gAlice->TreeR() ) {
349 sprintf(branchname,"%sPpsdRP",GetName()) ;
351 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
354 fTrackSegments->Delete() ;
356 if ( fTrackSegments && gAlice->TreeR() ) {
357 sprintf(branchname,"%sTS",GetName()) ;
358 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
361 fRecParticles->Delete() ;
363 if ( fRecParticles && gAlice->TreeR() ) {
364 sprintf(branchname,"%sRP",GetName()) ;
365 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
370 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
372 // 4. Expand or Shrink the arrays to the proper size
376 size = fEmcRecPoints->GetEntries() ;
377 fEmcRecPoints->Expand(size) ;
379 size = fPpsdRecPoints->GetEntries() ;
380 fPpsdRecPoints->Expand(size) ;
382 size = fTrackSegments->GetEntries() ;
383 fTrackSegments->Expand(size) ;
385 size = fRecParticles->GetEntries() ;
386 fRecParticles->Expand(size) ;
388 gAlice->TreeR()->Fill() ;
391 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
393 // Deleting reconstructed objects
394 ResetReconstruction();
398 //____________________________________________________________________________
399 void AliPHOSv1::ResetDigits()
401 // May sound strange, but cumulative hits are store in digits Tree
402 AliDetector::ResetDigits();
408 //____________________________________________________________________________
409 void AliPHOSv1::ResetReconstruction()
411 // Deleting reconstructed objects
413 if ( fEmcRecPoints ) fEmcRecPoints->Delete();
414 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
415 if ( fTrackSegments ) fTrackSegments->Delete();
416 if ( fRecParticles ) fRecParticles->Delete();
420 //____________________________________________________________________________
421 void AliPHOSv1::SetTreeAddress()
424 AliPHOS::SetTreeAddress();
426 // //Branch address for TreeR: RecPpsdRecPoint
427 // TTree *treeR = gAlice->TreeR();
428 // if ( treeR && fPpsdRecPoints ) {
429 // branch = treeR->GetBranch("PHOSPpsdRP");
430 // if (branch) branch->SetAddress(&fPpsdRecPoints) ;
434 //____________________________________________________________________________
436 void AliPHOSv1::StepManager(void)
438 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
440 Int_t relid[4] ; // (box, layer, row, column) indices
441 Float_t xyze[4] ; // position wrt MRS and energy deposited
445 Int_t tracknumber = gAlice->CurrentTrack() ;
446 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
447 TString name = fGeom->GetName() ;
448 if ( name == "GPS2" ) { // the CPV is a PPSD
449 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
451 gMC->TrackPosition(pos) ;
455 xyze[3] = gMC->Edep() ;
457 if ( xyze[3] != 0 ) { // there is deposited energy
458 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
459 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
460 // 1-> Geom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
461 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
462 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
463 gMC->CurrentVolID(relid[3]) ; // get the column number
465 // get the absolute Id number
468 fGeom->RelToAbsNumbering(relid, absid) ;
470 // add current hit to the hit list
471 AddHit(fIshunt, primary, tracknumber, absid, xyze);
473 } // there is deposited energy
474 } // We are inside the gas of the CPV
475 } // GPS2 configuration
477 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) // We are inside a PBWO crystal
479 gMC->TrackPosition(pos) ;
483 xyze[3] = gMC->Edep() ;
485 if ( xyze[3] != 0 ) {
486 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
487 relid[1] = 0 ; // means PBW04
488 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
489 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
491 // get the absolute Id number
494 fGeom->RelToAbsNumbering(relid, absid) ;
496 // add current hit to the hit list
498 AddHit(fIshunt, primary,tracknumber, absid, xyze);
500 } // there is deposited energy
501 } // we are inside a PHOS Xtal