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 **************************************************************************/
19 //_________________________________________________________________________
20 // Implementation version vImpacts of PHOS Manager class.
21 // This class inherits from v1 and adds impacts storing.
22 // Impacts stands for exact values of track coming to the detectors
24 // Impacts are written to the same tree as hits are
25 // but in separate branches.
27 //*-- Author: Yuri Kharlov (IHEP, Protvino/SUBATECH, Nantes)
30 // --- ROOT system ---
33 #include <TVirtualMC.h>
35 // --- Standard library ---
37 // --- AliRoot header files ---
39 #include "AliPHOSGeometry.h"
40 #include "AliPHOSImpact.h"
41 #include "AliPHOSvImpacts.h"
45 ClassImp(AliPHOSvImpacts)
47 //____________________________________________________________________________
48 AliPHOSvImpacts::AliPHOSvImpacts():AliPHOSv1()
56 //____________________________________________________________________________
57 AliPHOSvImpacts::AliPHOSvImpacts(const char *name, const char *title):
60 // ctor : title is used to identify the layout
63 // - fHits (the "normal" one), which retains the hits associated with
64 // the current primary particle being tracked
65 // (this array is reset after each primary has been tracked).
66 // This part inherits from AliPHOSv1
69 // - fEMCImpacts, fCPVImpacts which are
70 // TList of EMC and CPV modules respectively, each
71 // modules contains TClonesArray of AliPHOSImpacts
73 fEMCImpacts = new TList();
74 fCPVImpacts = new TList();
76 Int_t nPHOSModules = GetGeometry()->GetNModules();
77 Int_t nCPVModules = GetGeometry()->GetNModules();
80 TClonesArray * impacts;
81 for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
82 fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
83 fNEMCImpacts[iPHOSModule] = 0;
84 impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(iPHOSModule));
86 for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
87 fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
88 fNCPVImpacts[iPHOSModule] = 0;
89 impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(iPHOSModule));
94 //____________________________________________________________________________
95 AliPHOSvImpacts::~AliPHOSvImpacts()
106 // Delete impacts in EMC, CPV
108 fEMCImpacts->Delete() ;
113 fCPVImpacts->Delete() ;
119 //____________________________________________________________________________
120 void AliPHOSvImpacts::AddImpact( char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
121 Int_t pid, TLorentzVector p, Float_t *xyz)
123 // Add an impact to the impact list.
125 TClonesArray * impacts = 0;
128 if (strcmp(det,"EMC ")==0) {
129 impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(module));
130 nImpacts= fNEMCImpacts[module];
131 fNEMCImpacts[module]++ ;
133 else if (strcmp(det,"CPV ")==0) {
134 impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(module));
135 nImpacts= fNCPVImpacts[module];
136 fNCPVImpacts[module]++ ;
139 new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
142 printf("Module %d %s: ",module,det);
143 (dynamic_cast<AliPHOSImpact*>((impacts->At(nImpacts))))->Print();
147 //____________________________________________________________________________
148 void AliPHOSvImpacts::MakeBranch(Option_t *opt)
150 // Create new branch in the current Hits Root Tree containing
151 // a list of PHOS impacts (exact values of track coming to detector)
153 AliDetector::MakeBranch(opt);
155 Int_t bufferSize = 32000 ;
156 Int_t splitlevel = 0 ;
157 TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
158 TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
162 //____________________________________________________________________________
163 void AliPHOSvImpacts::ResetHits()
165 // Reset impact branches for EMC, CPV and PPSD
167 AliDetector::ResetHits();
170 for (i=0; i<GetGeometry()->GetNModules(); i++) {
171 (dynamic_cast<TClonesArray*>(fEMCImpacts->At(i))) -> Clear();
172 fNEMCImpacts[i] = 0 ;
175 for (i=0; i<GetGeometry()->GetNModules(); i++) {
176 (dynamic_cast<TClonesArray*>(fCPVImpacts->At(i))) -> Clear();
177 fNCPVImpacts[i] = 0 ;
182 //_____________________________________________________________________________
183 void AliPHOSvImpacts::StepManager(void)
185 // Find impacts (tracks which enter the EMC, CPV)
186 // and add them to to the impact lists
188 AliPHOSv1::StepManager();
190 Float_t xyzm[3], xyzd[3], pm[3], pd[3];
191 TLorentzVector pmom ; // Lorentz momentum of the particle initiated hit
192 TLorentzVector pos ; // Lorentz vector of the track current position
195 Int_t tracknumber = gAlice->GetMCApp()->GetCurrentTrackNumber() ;
196 Int_t primary = gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() );
197 TString name = GetGeometry()->GetName() ;
201 if( gMC->CurrentVolID(copy) == gMC->VolId("PXTL") &&
202 gMC->IsTrackEntering() ) {
203 gMC->TrackMomentum(pmom);
204 gMC->TrackPosition(pos) ;
207 for (i=0; i<3; i++) xyzm[i] = pos[i];
209 for (i=0; i<3; i++) {
213 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
214 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
216 // Select tracks coming to the crystal from up or down sides
217 if (pd[1]<0 && xyzd[1] > GetGeometry()->GetCrystalSize(1)/2-0.1 ||
218 pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.1) {
219 // Select tracks coming to the crystal from up or down sides
220 Int_t pid = gMC->TrackPid();
222 gMC->CurrentVolOffID(10,module);
224 AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
230 if( gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
231 gMC->IsTrackEntering() ) {
232 gMC->TrackMomentum(pmom);
233 gMC->TrackPosition(pos) ;
236 for (i=0; i<3; i++) xyzm[i] = pos[i];
238 for (i=0; i<3; i++) {
242 Int_t pid = gMC->TrackPid();
244 gMC->CurrentVolOffID(3,module);
246 AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);