]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSvImpacts.cxx
Create ideal, residual and fully displaced alignment objects in one macro
[u/mrichter/AliRoot.git] / PHOS / AliPHOSvImpacts.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 /* $Id$ */
18
19 /* History of cvs commits:
20  *
21  * $Log$
22  * Revision 1.21  2005/05/28 14:19:05  schutz
23  * Compilation warnings fixed by T.P.
24  *
25  */
26
27 //_________________________________________________________________________
28 // Implementation version vImpacts of PHOS Manager class.
29 // This class inherits from v1 and adds impacts storing.
30 // Impacts stands for exact values of track coming to the detectors
31 // EMC, CPV or PPSD.
32 // Impacts are written to the same tree as hits are
33 // but in separate branches.
34 //---
35 //*-- Author: Yuri Kharlov (IHEP, Protvino/SUBATECH, Nantes)
36
37
38 // --- ROOT system ---
39
40 //#include <TTree.h>
41 #include <TVirtualMC.h>
42
43 // --- Standard library ---
44
45 // --- AliRoot header files ---
46
47 #include "AliPHOSGeometry.h"
48 #include "AliPHOSImpact.h"
49 #include "AliPHOSvImpacts.h"
50 #include "AliRun.h"
51 #include "AliMC.h"
52 #include "AliLog.h"
53
54 ClassImp(AliPHOSvImpacts)
55
56 //____________________________________________________________________________
57 AliPHOSvImpacts::AliPHOSvImpacts():AliPHOSv1()
58 {
59   // ctor
60   fEMCImpacts  = 0 ;
61   fPPSDImpacts = 0 ;
62   fCPVImpacts  = 0 ; 
63 }
64
65 //____________________________________________________________________________
66 AliPHOSvImpacts::AliPHOSvImpacts(const char *name, const char *title):
67 AliPHOSv1(name,title) 
68 {
69   // ctor : title is used to identify the layout
70   //
71   // We store hits :
72   //   - fHits (the "normal" one), which retains the hits associated with
73   //     the current primary particle being tracked
74   //     (this array is reset after each primary has been tracked).
75   //     This part inherits from AliPHOSv1
76   //
77   // We store impacts :
78   //  - fEMCImpacts, fCPVImpacts which are
79   //    TList of EMC and CPV modules respectively, each
80   //    modules contains TClonesArray of AliPHOSImpacts
81   
82   fEMCImpacts  = new TList();
83   fCPVImpacts  = new TList();
84
85   Int_t nPHOSModules = GetGeometry()->GetNModules();
86   Int_t nCPVModules  = GetGeometry()->GetNModules();
87
88   Int_t iPHOSModule;
89   TClonesArray * impacts;
90   for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
91     fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
92     fNEMCImpacts[iPHOSModule] = 0;
93     impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(iPHOSModule));
94   }
95   for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
96     fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
97     fNCPVImpacts[iPHOSModule] = 0;
98     impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(iPHOSModule));
99   }
100
101 }
102
103 //____________________________________________________________________________
104 AliPHOSvImpacts::~AliPHOSvImpacts()
105 {
106   // dtor
107
108   // Delete hits
109   if ( fHits ) {
110     fHits->Delete() ; 
111     delete fHits ;
112     fHits = 0 ; 
113   }
114
115   // Delete impacts in EMC, CPV
116   if ( fEMCImpacts ) {
117     fEMCImpacts->Delete() ; 
118     delete fEMCImpacts ;
119     fEMCImpacts = 0 ; 
120   }
121   if ( fCPVImpacts ) {
122     fCPVImpacts->Delete() ; 
123     delete fCPVImpacts ;
124     fCPVImpacts = 0 ; 
125   }
126 }
127
128 //____________________________________________________________________________
129 void AliPHOSvImpacts::Copy(TObject & phos) const
130 {
131   TObject::Copy(phos) ; 
132   AliPHOS::Copy(phos) ; 
133 }
134
135
136 //____________________________________________________________________________
137 void AliPHOSvImpacts::AddImpact(const char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
138                            Int_t pid, TLorentzVector p, Float_t *xyz)
139 {
140   // Add an impact to the impact list.
141
142   TClonesArray * impacts = 0;
143   Int_t         nImpacts = 0;
144
145   if (strcmp(det,"EMC ")==0) {
146     impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(module));
147     nImpacts= fNEMCImpacts[module];
148     fNEMCImpacts[module]++ ;
149   }
150   else if (strcmp(det,"CPV ")==0) {
151     impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(module));
152     nImpacts= fNCPVImpacts[module];
153     fNCPVImpacts[module]++ ;
154   }
155
156   new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
157
158   AliDebugClass(1,Form("Module %d %s: ",module,det));
159   if (AliLog::GetGlobalDebugLevel()>0)
160     (dynamic_cast<AliPHOSImpact*>((impacts->At(nImpacts))))->Print();
161 }
162
163 //____________________________________________________________________________
164 void AliPHOSvImpacts::MakeBranch(Option_t *opt)
165 {  
166   // Create new branch in the current Hits Root Tree containing
167   // a list of PHOS impacts (exact values of track coming to detector)
168
169   AliDetector::MakeBranch(opt);
170   
171   Int_t bufferSize = 32000 ;
172   Int_t splitlevel = 0 ;
173   TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
174   TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
175   
176 }
177
178 //____________________________________________________________________________
179 void AliPHOSvImpacts::ResetHits() 
180 {              
181   // Reset impact branches for EMC, CPV and PPSD
182
183   AliDetector::ResetHits();
184
185   Int_t i;
186   for (i=0; i<GetGeometry()->GetNModules(); i++) {
187     (dynamic_cast<TClonesArray*>(fEMCImpacts->At(i))) -> Clear();
188     fNEMCImpacts[i] = 0 ;
189   }
190
191   for (i=0; i<GetGeometry()->GetNModules(); i++) {
192     (dynamic_cast<TClonesArray*>(fCPVImpacts->At(i))) -> Clear();
193     fNCPVImpacts[i] = 0 ;
194   }
195   
196 }
197
198 //_____________________________________________________________________________
199 void AliPHOSvImpacts::StepManager(void)
200 {
201   // Find impacts (tracks which enter the EMC, CPV)
202   // and add them to to the impact lists
203
204   AliPHOSv1::StepManager();
205
206   Float_t xyzm[3], xyzd[3], pm[3], pd[3];
207   TLorentzVector pmom     ;           // Lorentz momentum of the particle initiated hit
208   TLorentzVector pos      ;           // Lorentz vector of the track current position
209   Int_t          copy     ;
210
211   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber() ; 
212   Int_t primary     =  gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() ); 
213   TString name      =  GetGeometry()->GetName() ; 
214
215   // Add impact to EMC
216
217   static Int_t idPXTL = gMC->VolId("PXTL");
218   if( gMC->CurrentVolID(copy) == idPXTL &&
219       gMC->IsTrackEntering() ) {
220     gMC->TrackMomentum(pmom);
221     gMC->TrackPosition(pos) ;
222
223     Int_t i;
224     for (i=0; i<3; i++) xyzm[i] = pos[i];
225
226     for (i=0; i<3; i++) {
227       xyzm[i] = pos[i] ;
228       pm[i]   = pmom[i];
229     }
230     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
231     gMC -> Gmtod (pm,   pd,   2);    // transform 3-momentum from master to daughter system
232
233     // Select tracks coming to the crystal from up or down sides
234     if (pd[1]<0 && xyzd[1] >  GetGeometry()->GetCrystalSize(1)/2-0.1 ||
235         pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.1) {
236     // Select tracks coming to the crystal from up or down sides
237       Int_t pid = gMC->TrackPid();
238       Int_t module;
239       gMC->CurrentVolOffID(10,module);
240       module--;
241       AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
242     }
243   }
244
245   // Add impact to CPV
246
247   static Int_t idPCPQ = gMC->VolId("PCPQ");
248   if( gMC->CurrentVolID(copy) == idPCPQ &&
249       gMC->IsTrackEntering() ) {
250     gMC->TrackMomentum(pmom);
251     gMC->TrackPosition(pos) ;
252
253     Int_t i;
254     for (i=0; i<3; i++) xyzm[i] = pos[i];
255
256     for (i=0; i<3; i++) {
257       xyzm[i] = pos[i] ;
258       pm[i]   = pmom[i];
259     }
260     Int_t pid = gMC->TrackPid();
261     Int_t module;
262     gMC->CurrentVolOffID(3,module);
263     module--;
264     AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
265   }
266   
267 }