]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDv1.cxx
Fixed compilation warnings with gcc-4.
[u/mrichter/AliRoot.git] / FMD / AliFMDv1.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 /* $Id$ */
16 /** @file    AliFMDv1.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:48:51 2006
19     @brief   Concrete implementation of FMD detector driver - detailed
20     version 
21     @ingroup FMD_sim
22 */
23 //____________________________________________________________________
24 //                                                                          
25 // Forward Multiplicity Detector based on Silicon wafers. This class
26 // contains the base procedures for the Forward Multiplicity detector
27 // Detector consists of 3 sub-detectors FMD1, FMD2, and FMD3, each of
28 // which has 1 or 2 rings of silicon sensors. 
29 // This class contains the detailed version of the FMD - that is, hits
30 // are produced during simulation. 
31 //                                                                           
32 // See also the class AliFMD for a more detailed explanation of the
33 // various componets. 
34 //
35 #include <TVirtualMC.h>         // ROOT_TVirtualMC
36 #include <AliRun.h>             // ALIRUN_H
37 #include <AliMC.h>              // ALIMC_H
38 #include <AliLog.h>             // ALILOG_H
39 #include "AliFMDv1.h"           // ALIFMDV1_H
40 // #include "AliFMDGeometryBuilder.h"
41 #include "AliFMDGeometry.h"
42 #include "AliFMDDetector.h"
43 #include "AliFMDRing.h"
44 #include <TParticlePDG.h>
45 #include <TDatabasePDG.h>
46 #include "AliFMDHit.h"
47
48 //____________________________________________________________________
49 ClassImp(AliFMDv1)
50 #if 0
51   ; // This is here to keep Emacs for indenting the next line
52 #endif
53
54
55 //____________________________________________________________________
56 Bool_t
57 AliFMDv1::VMC2FMD(TLorentzVector& v, UShort_t& detector,
58                   Char_t& ring, UShort_t& sector, UShort_t& strip) const
59 {
60   // Convert VMC coordinates to detector coordinates 
61   TVirtualMC* mc = TVirtualMC::GetMC();
62   AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
63
64   // Get track position
65   mc->TrackPosition(v);
66   Int_t moduleno; mc->CurrentVolOffID(fmd->GetModuleOff(), moduleno);
67   Int_t iring;    mc->CurrentVolOffID(fmd->GetRingOff(), iring);   
68   ring = Char_t(iring);
69   Int_t det;      mc->CurrentVolOffID(fmd->GetDetectorOff(), det); 
70   detector = det;
71   
72
73   // Get the ring geometry
74   //Int_t     nsec = fmd->GetDetector(detector)->GetRing(ring)->GetNSectors();
75   Int_t     nstr  = fmd->GetDetector(detector)->GetRing(ring)->GetNStrips();
76   Double_t  lowr  = fmd->GetDetector(detector)->GetRing(ring)->GetMinR();
77   Double_t  theta = fmd->GetDetector(detector)->GetRing(ring)->GetTheta();
78   Double_t  pitch = fmd->GetDetector(detector)->GetRing(ring)->GetPitch();
79
80   // Figure out the strip number
81   Double_t r     = TMath::Sqrt(v.X() * v.X() + v.Y() * v.Y());
82   Int_t    str   = Int_t((r - lowr) / pitch);
83   if (str < 0 || str >= nstr) return kFALSE;
84   strip          = str;
85
86   // Figure out the sector number
87   Double_t phi    = TMath::ATan2(v.Y(), v.X()) * 180. / TMath::Pi();
88   if (phi < 0) phi = 360. + phi;
89   Double_t t      = phi - 2 * moduleno * theta;
90   sector          = 2 * moduleno;
91   if (t < 0 || t > 2 * theta) return kFALSE;
92   else if (t > theta)         sector += 1;
93
94   AliDebug(40, Form("<1> Inside an active FMD volume FMD%d%c[%2d,%3d] %s",
95                     detector, ring, sector, strip, mc->CurrentVolPath()));
96   return kTRUE;
97 }
98
99 //____________________________________________________________________
100 Bool_t
101 AliFMDv1::VMC2FMD(Int_t copy, TLorentzVector& v,
102                   UShort_t& detector, Char_t& ring,
103                   UShort_t& sector, UShort_t& strip) const
104 {
105   // Convert VMC coordinates to detector coordinates 
106   TVirtualMC* mc = TVirtualMC::GetMC();
107   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
108
109   strip = copy - 1;
110   Int_t sectordiv; mc->CurrentVolOffID(fmd->GetSectorOff(), sectordiv);
111   if (fmd->GetModuleOff() >= 0) {
112     Int_t module;    mc->CurrentVolOffID(fmd->GetModuleOff(), module);
113     sector = 2 * module + sectordiv;
114   }
115   else 
116     sector = sectordiv;
117   AliDebug(30, Form("Getting ring volume with offset %d -> %s", 
118                     fmd->GetRingOff(), 
119                     mc->CurrentVolOffName(fmd->GetRingOff())));
120   Int_t iring;     mc->CurrentVolOffID(fmd->GetRingOff(), iring); 
121   ring = Char_t(iring);
122   Int_t det;       mc->CurrentVolOffID(fmd->GetDetectorOff(), det); 
123   detector = det;
124
125   //Double_t  rz  = fmd->GetDetector(detector)->GetRingZ(ring);
126   AliFMDDetector* gdet  = fmd->GetDetector(detector);
127   AliFMDRing*     gring = gdet->GetRing(ring);
128   if (!gring) {
129     AliFatal(Form("Ring %c not found (volume was %s at offset %d in path %s)", 
130                   ring, fmd->GetRingOff(), 
131                   mc->CurrentVolOffName(fmd->GetRingOff()),
132                   mc->CurrentVolPath()));
133   }
134   Int_t n = gring->GetNSectors();
135 #if 0
136   if (rz < 0) {
137     Int_t s = ((n - sector + n / 2) % n) + 1;
138     AliDebug(1, Form("Recalculating sector to %d (=%d-%d+%d/2%%%d+1 z=%f)",
139                      s, n, sector, n, n, rz));
140     sector = s;
141   }
142 #endif
143   if (sector < 1 || sector > n) {
144     AliWarning(Form("sector # %d out of range (0-%d)", sector-1, n-1));
145     return kFALSE;
146   }
147   sector--;
148   // Get track position
149   mc->TrackPosition(v);
150   AliDebug(15, Form("<2> Inside an active FMD volume FMD%d%c[%2d,%3d] %s",
151                     detector, ring, sector, strip, mc->CurrentVolPath()));
152
153   return kTRUE;
154 }
155
156
157 //____________________________________________________________________
158 Bool_t
159 AliFMDv1::CheckHit(Int_t trackno, Int_t pdg, Float_t absQ, 
160                    const TLorentzVector& p, Float_t edep) const
161 {
162   // Check that a hit is good 
163   if (AliLog::GetDebugLevel("FMD", "AliFMD") < 5) return kFALSE;
164   TVirtualMC* mc   = TVirtualMC::GetMC();
165   Double_t mass    = mc->TrackMass();
166   Double_t poverm  = (mass == 0 ? 0 : p.P() / mass);
167     
168   // This `if' is to debug abnormal energy depositions.  We trigger on
169   // p/m approx larger than or equal to a MIP, and a large edep - more 
170   // than 1 keV - a MIP is 100 eV. 
171   if (!(edep > absQ * absQ && poverm > 1)) return kFALSE;
172   
173   TArrayI procs;
174   mc->StepProcesses(procs);
175   TString processes;
176   for (Int_t ip = 0; ip < procs.fN; ip++) {
177     if (ip != 0) processes.Append(",");
178     processes.Append(TMCProcessName[procs.fArray[ip]]);
179   }
180   TDatabasePDG* pdgDB        = TDatabasePDG::Instance();
181   TParticlePDG* particleType = pdgDB->GetParticle(pdg);
182   TString pname(particleType ? particleType->GetName() : "???");
183   TString what;
184   if (mc->IsTrackEntering())    what.Append("entering ");
185   if (mc->IsTrackExiting())     what.Append("exiting ");
186   if (mc->IsTrackInside())      what.Append("inside ");
187   if (mc->IsTrackDisappeared()) what.Append("disappeared ");
188   if (mc->IsTrackStop())        what.Append("stopped ");
189   if (mc->IsNewTrack())         what.Append("new ");
190   if (mc->IsTrackAlive())       what.Append("alive ");
191   if (mc->IsTrackOut())         what.Append("out ");
192       
193   Int_t mother = gAlice->GetMCApp()->GetPrimary(trackno);
194   AliDebug(15, Form("Track # %5d deposits a lot of energy\n" 
195                     "  Volume:    %s\n" 
196                     "  Momentum:  (%7.4f,%7.4f,%7.4f)\n"
197                     "  PDG:       %d (%s)\n" 
198                     "  Edep:      %-14.7f keV (mother %d)\n"
199                     "  p/m:       %-7.4f/%-7.4f = %-14.7f\n"
200                     "  Processes: %s\n"
201                     "  What:      %s\n",
202                     trackno, mc->CurrentVolPath(), p.X(), p.Y(), p.Z(),
203                     pdg, pname.Data(), edep, mother, p.P(), mass, 
204                     poverm, processes.Data(), what.Data()));
205   return kTRUE;
206 }
207
208
209 //____________________________________________________________________
210 void 
211 AliFMDv1::StepManager()
212 {
213   // Member function that is executed each time a hit is made in the
214   // FMD.  None-charged particles are ignored.   Dead tracks  are
215   // ignored. 
216   //
217   // The procedure is as follows: 
218   // 
219   //   - IF NOT track is alive THEN RETURN ENDIF
220   //   - IF NOT particle is charged THEN RETURN ENDIF
221   //   - IF NOT volume name is "STRI" or "STRO" THEN RETURN ENDIF 
222   //   - Get strip number (volume copy # minus 1)
223   //   - Get phi division number (mother volume copy #)
224   //   - Get module number (grand-mother volume copy #)
225   //   - section # = 2 * module # + phi division # - 1
226   //   - Get ring Id from volume name 
227   //   - Get detector # from grand-grand-grand-mother volume name 
228   //   - Get pointer to sub-detector object. 
229   //   - Get track position 
230   //   - IF track is entering volume AND track is inside real shape THEN
231   //   -   Reset energy deposited 
232   //   -   Get track momentum 
233   //   -   Get particle ID # 
234   ///  - ENDIF
235   //   - IF track is inside volume AND inside real shape THEN 
236   ///  -   Update energy deposited 
237   //   - ENDIF 
238   //   - IF track is inside real shape AND (track is leaving volume,
239   //         or it died, or it is stopped  THEN
240   //   -   Create a hit 
241   //   - ENDIF
242   //     
243   TVirtualMC* mc = TVirtualMC::GetMC();
244   if (!mc->IsTrackAlive()) return;
245   Double_t absQ = TMath::Abs(mc->TrackCharge());
246   if (absQ <= 0) return;
247   
248   Int_t copy;
249   Int_t vol = mc->CurrentVolID(copy);
250   AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
251   if (!fmd->IsActive(vol)) {
252     AliDebug(50, Form("Not an FMD volume %d '%s'",vol,mc->CurrentVolName()));
253     return;
254   }
255   TLorentzVector v;
256   UShort_t       detector;
257   Char_t         ring;
258   UShort_t       sector;
259   UShort_t       strip;
260   
261   if (fmd->IsDetailed()) {
262     if (!VMC2FMD(copy, v, detector, ring, sector, strip)) return;
263   } else {
264     if (!VMC2FMD(v, detector, ring, sector, strip)) return;
265   }
266   TLorentzVector p;
267   mc->TrackMomentum(p);
268   Int_t    trackno = gAlice->GetMCApp()->GetCurrentTrackNumber();
269   Int_t    pdg     = mc->TrackPid();
270   Double_t edep    = mc->Edep() * 1000; // keV
271   Bool_t   isBad   = CheckHit(trackno, pdg, absQ, p, edep);
272   
273   // Check that the track is actually within the active area 
274   Bool_t entering = mc->IsTrackEntering();
275   Bool_t inside   = mc->IsTrackInside();
276   Bool_t out      = (mc->IsTrackExiting()|| mc->IsTrackDisappeared()||
277                      mc->IsTrackStop());
278   // Reset the energy deposition for this track, and update some of
279   // our parameters.
280   if (entering) {
281     AliDebug(15, Form("Track # %8d entering active FMD volume %s: "
282                       "Edep=%f (%f,%f,%f)", trackno, mc->CurrentVolPath(),
283                       edep, v.X(), v.Y(), v.Z()));
284     fCurrentP      = p;
285     fCurrentV      = v;    
286     fCurrentDeltaE = edep;
287     fCurrentPdg    = pdg; // mc->IdFromPDG(pdg);
288   }
289   // If the track is inside, then update the energy deposition
290   if (inside && fCurrentDeltaE >= 0) {
291     fCurrentDeltaE += edep;
292     AliDebug(15, Form("Track # %8d inside active FMD volume %s: Edep=%f, "
293                       "Accumulated Edep=%f  (%f,%f,%f)", trackno, 
294                       mc->CurrentVolPath(), edep, fCurrentDeltaE, 
295                       v.X(), v.Y(), v.Z()));
296   }
297   // The track exits the volume, or it disappeared in the volume, or
298   // the track is stopped because it no longer fulfills the cuts
299   // defined, then we create a hit. 
300   if (out) {
301     if (fCurrentDeltaE >= 0) {
302       fCurrentDeltaE += edep;
303       AliDebug(15, Form("Track # %8d exiting active FMD volume %s: Edep=%g, "
304                         "Accumulated Edep=%g (%f,%f,%f)", trackno, 
305                         mc->CurrentVolPath(), edep, fCurrentDeltaE, 
306                         v.X(), v.Y(), v.Z()));
307       TVector3 cur(v.Vect());
308       cur -= fCurrentV.Vect();
309       Double_t len = cur.Mag();
310       AliFMDHit* h = 
311         AddHitByFields(trackno, detector, ring, sector, strip,
312                        fCurrentV.X(),  fCurrentV.Y(), fCurrentV.Z(),
313                        fCurrentP.X(),  fCurrentP.Y(), fCurrentP.Z(), 
314                        fCurrentDeltaE, fCurrentPdg,   fCurrentV.T(), 
315                        len, mc->IsTrackDisappeared()||mc->IsTrackStop());
316       // Add a copy 
317       if (isBad && fBad) { 
318         new ((*fBad)[fBad->GetEntries()]) AliFMDHit(*h);
319       }
320       // Check the geometry that we can get back the coordinates. 
321 #ifdef CHECK_TRANS
322       Double_t x, y, z;
323       fmd->Detector2XYZ(detector, ring, sector, strip, x, y ,z);
324       AliDebug(1, Form("Hit at (%f,%f,%f), geometry says (%f,%f,%f)", 
325                        fCurrentV.X(), fCurrentV.Y(), fCurrentV.Z(), x, y, z));
326 #endif
327     }
328     fCurrentDeltaE = -1;
329   }
330 }
331 //___________________________________________________________________
332 //
333 // EOF
334 //