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 vImpacts of PHOS Manager class.
20 // This class inherits from v1 and adds impacts storing.
21 // Impacts stands for exact values of track coming to the detectors
23 // Impacts are written to the same tree as hits are
24 // but in separate branches.
26 //*-- Author: Yuri Kharlov (IHEP, Protvino/SUBATECH, Nantes)
29 // --- ROOT system ---
33 // --- Standard library ---
35 // --- AliRoot header files ---
39 #include "AliPHOSvImpacts.h"
40 #include "AliPHOSGeometry.h"
41 #include "AliPHOSImpact.h"
43 ClassImp(AliPHOSvImpacts)
45 //____________________________________________________________________________
46 AliPHOSvImpacts::AliPHOSvImpacts()
51 //____________________________________________________________________________
52 AliPHOSvImpacts::AliPHOSvImpacts(const char *name, const char *title):
55 // ctor : title is used to identify the layout
56 // GPS2 = 5 modules (EMC + PPSD)
57 // IHEP = 5 modules (EMC + CPV )
58 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
61 // - fHits (the "normal" one), which retains the hits associated with
62 // the current primary particle being tracked
63 // (this array is reset after each primary has been tracked).
64 // This part inherits from AliPHOSv1
67 // - fEMCImpacts, fCPVImpacts, fPPSDImpacts which are
68 // TList of EMC, CPV and PPSD modules respectively, each
69 // modules contains TClonesArray of AliPHOSImpacts
71 fEMCImpacts = new TList();
72 fCPVImpacts = new TList();
73 fPPSDImpacts = new TList();
75 Int_t nPHOSModules = GetGeometry()->GetNModules();
76 Int_t nCPVModules = GetGeometry()->GetNCPVModules();
77 Int_t nPPSDModules = GetGeometry()->GetNPPSDModules();
80 TClonesArray * impacts;
81 for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
82 fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
83 fNEMCImpacts[iPHOSModule] = 0;
84 impacts = (TClonesArray *)fEMCImpacts->At(iPHOSModule);
86 for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
87 fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
88 fNCPVImpacts[iPHOSModule] = 0;
89 impacts = (TClonesArray *)fCPVImpacts->At(iPHOSModule);
91 for (iPHOSModule=0; iPHOSModule<nPPSDModules; iPHOSModule++) {
92 fPPSDImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
93 fNPPSDImpacts[iPHOSModule] = 0;
94 impacts = (TClonesArray *)fPPSDImpacts->At(iPHOSModule);
99 //____________________________________________________________________________
100 AliPHOSvImpacts::~AliPHOSvImpacts()
111 // Delete impacts in EMC, CPV and PPSD
113 fEMCImpacts->Delete() ;
118 fCPVImpacts->Delete() ;
122 if ( fPPSDImpacts ) {
123 fPPSDImpacts->Delete() ;
124 delete fPPSDImpacts ;
129 //____________________________________________________________________________
130 void AliPHOSvImpacts::AddImpact( char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
131 Int_t pid, TLorentzVector p, Float_t *xyz)
133 // Add an impact to the impact list.
135 TClonesArray * impacts = 0;
138 if (strcmp(det,"EMC ")==0) {
139 impacts = (TClonesArray *)fEMCImpacts->At(module);
140 nImpacts= fNEMCImpacts[module];
141 fNEMCImpacts[module]++ ;
143 else if (strcmp(det,"CPV ")==0) {
144 impacts = (TClonesArray *)fCPVImpacts->At(module);
145 nImpacts= fNCPVImpacts[module];
146 fNCPVImpacts[module]++ ;
148 else if (strcmp(det,"PPSD")==0) {
149 impacts = (TClonesArray *)fPPSDImpacts->At(module);
150 nImpacts= fNPPSDImpacts[module];
151 fNPPSDImpacts[module]++ ;
154 new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
157 printf("Module %d %s: ",module,det);
158 ((AliPHOSImpact*)(impacts->At(nImpacts)))->Print();
162 //____________________________________________________________________________
163 void AliPHOSvImpacts::MakeBranch(Option_t *opt, const char *file)
165 // Create new branch in the current Hits Root Tree containing
166 // a list of PHOS impacts (exact values of track coming to detector)
168 AliDetector::MakeBranch(opt,file);
170 Int_t bufferSize = 32000 ;
171 Int_t splitlevel = 0 ;
172 gAlice->TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
173 gAlice->TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
174 gAlice->TreeH()->Branch("PHOSPpsdImpacts", "TList", &fPPSDImpacts, bufferSize, splitlevel);
178 //____________________________________________________________________________
179 void AliPHOSvImpacts::ResetHits()
181 // Reset impact branches for EMC, CPV and PPSD
183 AliDetector::ResetHits();
186 for (i=0; i<GetGeometry()->GetNModules(); i++) {
187 ((TClonesArray*)fEMCImpacts->At(i)) -> Clear();
188 fNEMCImpacts[i] = 0 ;
191 if ( strcmp(GetGeometry()->GetName(),"IHEP") == 0 || strcmp(GetGeometry()->GetName(),"MIXT") == 0 ) {
192 for (i=0; i<GetGeometry()->GetNCPVModules(); i++) {
193 ((TClonesArray*)fCPVImpacts->At(i)) -> Clear();
194 fNCPVImpacts[i] = 0 ;
198 if ( strcmp(GetGeometry()->GetName(),"GPS2") == 0 || strcmp(GetGeometry()->GetName(),"MIXT") == 0 ) {
199 for (i=0; i<GetGeometry()->GetNPPSDModules(); i++) {
200 ((TClonesArray*)fPPSDImpacts->At(i)) -> Clear();
201 fNPPSDImpacts[i] = 0 ;
207 //_____________________________________________________________________________
208 void AliPHOSvImpacts::StepManager(void)
210 // Find impacts (tracks which enter the EMC, CPV or PPSD)
211 // and add them to to the impact lists
213 AliPHOSv1::StepManager();
215 Float_t xyzm[3], xyzd[3], pm[3], pd[3];
216 TLorentzVector pmom ; // Lorentz momentum of the particle initiated hit
217 TLorentzVector pos ; // Lorentz vector of the track current position
220 Int_t tracknumber = gAlice->CurrentTrack() ;
221 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
222 TString name = GetGeometry()->GetName() ;
226 if( gMC->CurrentVolID(copy) == gMC->VolId("PXTL") &&
227 gMC->IsTrackEntering() ) {
228 gMC->TrackMomentum(pmom);
229 gMC->TrackPosition(pos) ;
232 for (i=0; i<3; i++) xyzm[i] = pos[i];
234 for (i=0; i<3; i++) {
238 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
239 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
241 // Select tracks coming to the crystal from up or down sides
242 if (pd[1]<0 && xyzd[1] > GetGeometry()->GetCrystalSize(1)/2-0.001 ||
243 pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.001) {
244 Int_t pid = gMC->TrackPid();
246 gMC->CurrentVolOffID(10,module);
247 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
248 module += GetGeometry()->GetNModules() - GetGeometry()->GetNPPSDModules();
250 AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
256 if( (name == "IHEP" || name == "MIXT") &&
257 gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
258 gMC->IsTrackEntering() ) {
259 gMC->TrackMomentum(pmom);
260 gMC->TrackPosition(pos) ;
263 for (i=0; i<3; i++) xyzm[i] = pos[i];
265 for (i=0; i<3; i++) {
269 Int_t pid = gMC->TrackPid();
271 gMC->CurrentVolOffID(3,module);
273 AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
276 // Add impact to PPSD
278 if( (name == "GPS2" || name == "MIXT") &&
279 gMC->CurrentVolID(copy) == gMC->VolId("PPCE") &&
280 gMC->IsTrackEntering() ) {
281 gMC->TrackMomentum(pmom);
282 gMC->TrackPosition(pos) ;
285 for (i=0; i<3; i++) xyzm[i] = pos[i];
287 for (i=0; i<3; i++) {
291 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
292 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
294 // Select tracks coming to the crystal from up or down sides
295 if (pd[1]<0 && xyzd[1] > (GetGeometry()->GetConversionGap() + GetGeometry()->GetAvalancheGap())/2-0.001 ||
296 pd[1]>0 && xyzd[1] < -(GetGeometry()->GetConversionGap() + GetGeometry()->GetAvalancheGap())/2+0.001) {
297 Int_t pid = gMC->TrackPid();
299 gMC->CurrentVolOffID(5,module);
301 AddImpact("PPSD",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);