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