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