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 /* History of cvs commits:
22 * Revision 1.23 2006/09/13 07:31:01 kharlov
23 * Effective C++ corrections (T.Pocheptsov)
25 * Revision 1.22 2005/06/17 07:39:07 hristov
26 * Removing GetDebug and SetDebug from AliRun and AliModule. Using AliLog for the messages
28 * Revision 1.21 2005/05/28 14:19:05 schutz
29 * Compilation warnings fixed by T.P.
33 //_________________________________________________________________________
34 // Implementation version vImpacts of PHOS Manager class.
35 // This class inherits from v1 and adds impacts storing.
36 // Impacts stands for exact values of track coming to the detectors
38 // Impacts are written to the same tree as hits are
39 // but in separate branches.
41 //*-- Author: Yuri Kharlov (IHEP, Protvino/SUBATECH, Nantes)
44 // --- ROOT system ---
47 #include <TVirtualMC.h>
49 // --- Standard library ---
51 // --- AliRoot header files ---
53 #include "AliPHOSGeometry.h"
54 #include "AliPHOSImpact.h"
55 #include "AliPHOSvImpacts.h"
60 ClassImp(AliPHOSvImpacts)
62 //____________________________________________________________________________
63 AliPHOSvImpacts::AliPHOSvImpacts():
71 //____________________________________________________________________________
72 AliPHOSvImpacts::AliPHOSvImpacts(const char *name, const char *title):
73 AliPHOSv1(name,title),
74 fEMCImpacts(new TList),
75 fCPVImpacts(new TList),
78 // ctor : title is used to identify the layout
81 // - fHits (the "normal" one), which retains the hits associated with
82 // the current primary particle being tracked
83 // (this array is reset after each primary has been tracked).
84 // This part inherits from AliPHOSv1
87 // - fEMCImpacts, fCPVImpacts which are
88 // TList of EMC and CPV modules respectively, each
89 // modules contains TClonesArray of AliPHOSImpacts
91 Int_t nPHOSModules = GetGeometry()->GetNModules();
92 Int_t nCPVModules = GetGeometry()->GetNModules();
95 TClonesArray * impacts;
96 for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
97 fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
98 fNEMCImpacts[iPHOSModule] = 0;
99 impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(iPHOSModule));
101 for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
102 fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
103 fNCPVImpacts[iPHOSModule] = 0;
104 impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(iPHOSModule));
109 //____________________________________________________________________________
110 AliPHOSvImpacts::~AliPHOSvImpacts()
121 // Delete impacts in EMC, CPV
123 fEMCImpacts->Delete() ;
128 fCPVImpacts->Delete() ;
134 //____________________________________________________________________________
135 void AliPHOSvImpacts::AddImpact(const char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
136 Int_t pid, TLorentzVector p, Float_t *xyz)
138 // Add an impact to the impact list.
140 TClonesArray * impacts = 0;
143 if (strcmp(det,"EMC ")==0) {
144 impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(module));
145 nImpacts= fNEMCImpacts[module];
146 fNEMCImpacts[module]++ ;
148 else if (strcmp(det,"CPV ")==0) {
149 impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(module));
150 nImpacts= fNCPVImpacts[module];
151 fNCPVImpacts[module]++ ;
154 new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
156 AliDebugClass(1,Form("Module %d %s: ",module,det));
157 if (AliLog::GetGlobalDebugLevel()>0)
158 (dynamic_cast<AliPHOSImpact*>((impacts->At(nImpacts))))->Print();
161 //____________________________________________________________________________
162 void AliPHOSvImpacts::MakeBranch(Option_t *opt)
164 // Create new branch in the current Hits Root Tree containing
165 // a list of PHOS impacts (exact values of track coming to detector)
167 AliDetector::MakeBranch(opt);
169 Int_t bufferSize = 32000 ;
170 Int_t splitlevel = 0 ;
171 TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
172 TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
176 //____________________________________________________________________________
177 void AliPHOSvImpacts::ResetHits()
179 // Reset impact branches for EMC, CPV and PPSD
181 AliDetector::ResetHits();
184 for (i=0; i<GetGeometry()->GetNModules(); i++) {
185 (dynamic_cast<TClonesArray*>(fEMCImpacts->At(i))) -> Clear();
186 fNEMCImpacts[i] = 0 ;
189 for (i=0; i<GetGeometry()->GetNModules(); i++) {
190 (dynamic_cast<TClonesArray*>(fCPVImpacts->At(i))) -> Clear();
191 fNCPVImpacts[i] = 0 ;
196 //_____________________________________________________________________________
197 void AliPHOSvImpacts::StepManager(void)
199 // Find impacts (tracks which enter the EMC, CPV)
200 // and add them to to the impact lists
202 AliPHOSv1::StepManager();
204 Float_t xyzm[3], xyzd[3], pm[3], pd[3];
205 TLorentzVector pmom ; // Lorentz momentum of the particle initiated hit
206 TLorentzVector pos ; // Lorentz vector of the track current position
209 Int_t tracknumber = gAlice->GetMCApp()->GetCurrentTrackNumber() ;
210 Int_t primary = gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() );
211 TString name = GetGeometry()->GetName() ;
215 static Int_t idPXTL = gMC->VolId("PXTL");
216 if( gMC->CurrentVolID(copy) == idPXTL &&
217 gMC->IsTrackEntering() ) {
218 gMC->TrackMomentum(pmom);
219 gMC->TrackPosition(pos) ;
222 for (i=0; i<3; i++) xyzm[i] = pos[i];
224 for (i=0; i<3; i++) {
228 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
229 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
231 // Select tracks coming to the crystal from up or down sides
232 if (pd[1]<0 && xyzd[1] > GetGeometry()->GetCrystalSize(1)/2-0.1 ||
233 pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.1) {
234 // Select tracks coming to the crystal from up or down sides
235 Int_t pid = gMC->TrackPid();
237 gMC->CurrentVolOffID(10,module);
239 AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
245 static Int_t idPCPQ = gMC->VolId("PCPQ");
246 if( gMC->CurrentVolID(copy) == idPCPQ &&
247 gMC->IsTrackEntering() ) {
248 gMC->TrackMomentum(pmom);
249 gMC->TrackPosition(pos) ;
252 for (i=0; i<3; i++) xyzm[i] = pos[i];
254 for (i=0; i<3; i++) {
258 Int_t pid = gMC->TrackPid();
260 gMC->CurrentVolOffID(3,module);
262 AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);