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