All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
VTXDigitizer.h
Go to the documentation of this file.
1 /* -*- Mode: C++; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 #ifndef VTXDigitizer_h
3 #define VTXDigitizer_h 1
4 
5 #include "marlin/Processor.h"
6 #include "lcio.h"
7 #include "EVENT/SimTrackerHit.h"
8 #include "IMPL/TrackerHitImpl.h"
9 #include "IMPL/SimTrackerHitImpl.h"
10 #include <string>
11 #include <vector>
13 #include "EVENT/LCIO.h"
14 #include <IMPL/LCCollectionVec.h>
15 
16 
17 using namespace lcio ;
18 using namespace marlin ;
19 
20 /**
21 \addtogroup TrackDigi TrackDigi
22 @{
23 */
24 struct IonisationPoint {
25  double x{};
26  double y{};
27  double z{};
28  double eloss{};
29 };
30 
31 struct SignalPoint {
32  double x{};
33  double y{};
34  double sigmaX{};
35  double sigmaY{};
36  double charge{};
37 
38 };
39 
40 typedef std::vector<IonisationPoint> IonisationPointVec;
41 typedef std::vector<SignalPoint> SignalPointVec;
42 typedef std::vector<TrackerHitImpl*> TrackerHitImplVec;
43 typedef std::vector<SimTrackerHitImpl*> SimTrackerHitImplVec;
44 
45 /**
46 \addtogroup VTXDigitizer VTXDigitizer
47 @{
48 Digitizer for Simulated Hits in the Vertex Detector. <br>
49  * Digitization follows the procedure adopted in the CMS software package. <br>
50  * For each Simulated Tracker Hit the intersection points with <br>
51  * the inner and outer boundaries of sensitive Si layer are calculated. <br>
52  * Track segment within a layer is approximated with the line. <br>
53  * It is divided by n subsegments, where n can be specified by a user <br>
54  * via external Processor parameter. For each subsegment the charge <br>
55  * is simulated according to Landau distribution as expected for Silicon. <br>
56  * The charge transfer from the middle point of track subsegment (referred <br>
57  * hereafter to as ionisation point) <br>
58  * to the outer collection plane is performed taking <br>
59  * into account Lorentz effect in the magnetic field <br>
60  * It is assumed that on the collection plane the electron cloud from <br>
61  * each ionisation point is spread according to the Gaussian distribution <br>
62  * whose width is proportional to the square-root of the drift distance. <br>
63  * The charge on each fired pixel is then calculated as a sum of contributions <br>
64  * from n Gaussians. The VTX ladder is assumed to have matrix of rectangular pixels <br>
65  * (In the future I plan to implement possibility of variying pixel size with <br>
66  * the z coordinate). <br>
67  * The output of the processor is the collection of Reconstructed Tracker Hits. <br>
68  * Each reconstructed hit is assigned the position of the center-of-gravity of <br>
69  * the cluster of fired pixels. Lorentz effect is corrected for.
70  * <h4>Input collections and prerequisites</h4>
71  * Processor requires collection of simulated vertex tracker hits. <br>
72  * If such a collection with the user specified name does not exist <br>
73  * processor takes no action.
74  * <h4>Output</h4>
75  * Processor produces an output collection of the Tracker Hits.
76  * Collection has name "VTXTrackerHits".
77  * @param CollectionName name of input SimTrackerHit collection <br>
78  * (default parameter value : "vxd01_VXD", taken from Mokka)
79  * @param TanLorentz tangent of the Lorentz angle <br>
80  * (default parameter value : 0.8) <br>
81  * @param CutOnDeltaRays cut on the energy of delta-electrons (in MeV) <br>
82  * (default parameter value : 0.03) <br>
83  * @param Diffusion diffusion coefficient for the nominal active layer thickness (in mm) <br>
84  * (default parameter value : 0.002) <br>
85  * @param LayerThickness thickness of the active Silicon layer (in mm) <br>
86  * (default parameter value : 0.03744) <br>
87  * @param PixelSizeX pixel size along direction perpendicular to beam axis (in mm) <br>
88  * (default value : 0.025) <br>
89  * @param PixelSizeY pixel size along beam axis (in mm) <br>
90  * (default value : 0.025) <br>
91  * @param ElectronsPerMeV number of electrons produced per MeV of deposited energy <br>
92  * (default parameter value : 270.3) <br>
93  * @param Threshold threshold on charge deposited on one pixel (in electons) <br>
94  * (default parameter value : 200.0) <br>
95  * @param LaddersInLayer vector of integers, numbers of phi-ladders in each layer <br>
96  * (default parameter values : 8, 8, 12, 16, 20; taken from Mokka database for VXD model vxd01) <br>
97  * @param LadderRadius vector of doubles, radii of layers (in mm) <br>
98  * (default parameter values : 15.301, 26.301, 38.301, 49.301, 60.301; taken from Mokka database
99  * for VXD model vxd01) <br>
100  * @param ActiveLadderOffset vector of doubles, offset of active Si layer along phi-angle (in mm) for each layer <br>
101  * (default parameter values : 1.455, 1.39866, 2.57163, 3.59295, 4.42245) <br>
102  * @param LadderHalfWidth vector of doubles, half-width of the ladders in each layer (in mm)<br>
103  * (default parameter values : 6.5, 11.0, 11.0, 11.0, 11.0; taken from Mokka database for VXD model vxd01) <br>
104  * @param LadderGaps vector of doubles, gaps between two symmetric ladders (+/-z) in mm<br>
105  * (default parameter values : 0.0, 0.04, 0.04, 0.04, 0.04; taken from Mokka database) <br>
106  * @param PhiOffset vector of doubles, offset in phi angle for starting ladder in each layer <br>
107  * (default parameter values : 0, 0, 0, 0, 0; taken from Mokka database for VXD model vxd01) <br>
108  * @param SegmentLength segment length along track path which is used to subdivide track into segments (in mm).
109  * The number of track subsegments is calculated as int(TrackLengthWithinActiveLayer/SegmentLength)+1 <br>
110  * (default parameter value : 0.005) <br>
111  * @param WidthOfCluster defines width in Gaussian sigmas to perform charge integration for
112  * a given pixel <br>
113  * (default parameter value : 3.0) <br>
114  * @param ElectronicEffects flag to switch on gaussian smearing of signal (electronic noise) <br>
115  * (default parameter value : 1) <br>
116  * @param ElectronicNoise electronic noise in electrons <br>
117  * (default parameter value : 100) <br>
118  * @param GenerateBackground flag to switch on pseudo-generation of pair background hits
119  * Background hits are uniformly generated in cosQ and Phi <br>
120  * (default parameter value : 0)
121  * @param BackgroundHitsPerLayer expected mean value of background hits accumulated in each layer
122  * over readout time. This number is calculated as Number of background hits per bunch crossing times
123  * number of bunch crossings over integration time. <br>
124  * (default values : 34400 23900 9600 5500 3100 corresponding to 100 overlaid bunch crossings) <br>
125  * @param Debug boolean variable, if set to 1, printout is activated <br>
126  * (default parameter value : 0) <br>
127  * <br>
128  * @author A. Raspereza, MPI Munich
129  */
130 class VTXDigitizer : public Processor {
131 
132  public:
133  VTXDigitizer(const VTXDigitizer&) = delete;
134  VTXDigitizer& operator=(const VTXDigitizer&) = delete;
135 
136  virtual Processor* newProcessor() { return new VTXDigitizer ; }
137 
138 
139  VTXDigitizer() ;
140 
141  /**
142  * Initialisation member function
143  */
144  virtual void init() ;
145 
146  /**
147  * Processing of run header
148  */
149  virtual void processRunHeader( LCRunHeader* run ) ;
150 
151  /** Processing of one event
152  */
153  virtual void processEvent( LCEvent * evt ) ;
154 
155  /** Produces check plots
156  */
157  virtual void check( LCEvent * evt ) ;
158 
159  /** Termination member function
160  */
161  virtual void end() ;
162 
163 
164  protected:
165 
166  /** Input collection name.
167  */
168  std::string _colName{};
169  std::string _outputCollectionName{};
170  std::string _colVTXRelation{};
171 
172  /** Run number
173  */
174  int _nRun{};
175 
176  /** Event number
177  */
178  int _nEvt{};
179  /** tangent of Lorentz angle
180  */
181  double _tanLorentzAngle{};
182  /** cut in MeV on delta electrons
183  * used in simulation of charge for each ionisation point
184  */
185  double _cutOnDeltaRays{};
186  /** Diffusion coefficient in mm for nominla layer thickness
187  */
188  double _diffusionCoefficient{};
189  /** layer thickness
190  */
191  // double _layerThickness;
192 // /** layer half-thickness
193 // */
194 // double _layerHalfThickness;
195 
196 
197  int _numberOfLayers{};
198  double _pixelSizeX{};
199  double _pixelSizeY{};
200  double _electronsPerKeV{};
201  double _segmentDepth{};
202  double _currentTotalCharge{};
203 
204  std::vector<int> _laddersInLayer{};
205  std::vector<float> _layerRadius{};
206  std::vector<float> _layerThickness{};
207  std::vector<float> _layerHalfThickness{};
208  std::vector<float> _layerLadderLength{};
209  std::vector<float> _layerLadderHalfWidth{};
210  std::vector<float> _layerPhiOffset{};
211  std::vector<float> _layerActiveSiOffset{};
212  std::vector<float> _layerHalfPhi{};
213  std::vector<float> _layerLadderGap{};
214  std::vector<float> _bkgdHitsInLayer{};
215  std::vector<float> _layerLadderWidth{};
216 
217  int _currentLayer{};
218  int _currentModule{};
219  int _generateBackground{};
220  double _currentParticleMomentum{};
221  double _currentParticleEnergy{};
222  double _currentParticleMass{};
223  double _currentPhi{};
224  double _widthOfCluster{};
225 
226  double PI{},TWOPI{},PI2{};
227  double SCALING{};
228 
229  int _produceFullPattern{};
230  int _numberOfSegments{};
231  int _debug{};
232  int _PoissonSmearing{};
233  int _electronicEffects{};
234  int _useMCPMomentum{};
235  int _removeDrays{};
236 
237  double _threshold{};
238  double _currentLocalPosition[3]{};
239  double _currentEntryPoint[3]{};
240  double _currentExitPoint[3]{};
241  double _electronicNoise{};
242  double _segmentLength{};
243 
244  IonisationPointVec _ionisationPoints{};
245  SignalPointVec _signalPoints{};
246 
248 
249  /**
250  * Finds coordinates
251  */
252  void FindLocalPosition(SimTrackerHit * hit,
253  double * localPosition,
254  double * localDirection);
255 
256  void TransformToLab(double * xLoc, double * xLab);
257  void ProduceIonisationPoints( SimTrackerHit * hit);
258  void ProduceSignalPoints( );
259  void ProduceHits(SimTrackerHitImplVec & simTrkVec);
260  void TransformXYToCellID(double x, double y,
261  int & ix,
262  int & iy);
263  void TransformCellIDToXY(int ix, int iy,
264  double & x, double & y);
265  void PoissonSmearer( SimTrackerHitImplVec & simTrkVec );
266  void GainSmearer( SimTrackerHitImplVec & simTrkVec );
267  void PrintInfo( SimTrackerHit * simTrkHit, TrackerHitImpl * recoHit);
268  TrackerHitImpl * ReconstructTrackerHit(SimTrackerHitImplVec & simTrkVec );
269  void TrackerHitToLab( TrackerHitImpl * recoHit );
270  void PositionWithinCell(double x, double y,
271  int & ix, int & iy,
272  double & xCell, double & yCell);
273  void generateBackground(LCCollectionVec * col);
274 
275  double _xLayerReco{},_yLayerReco{},_zLayerReco{};
276  double _xLayerSim{},_yLayerSim{},_zLayerSim{};
277  int _iLayer{};
278 
279  int _nCoveredX{},_nCoveredY{},_nCells{};
280 
281  int _totEntries{};
282  double _totMomentum{};
283  double _momX{},_momY{},_momZ{};
284  double _eDep{};
285 
286  double _amplX[20]{};
287  double _amplY[20]{};
288  double _amplC[100]{};
289  double _ampl{};
290  double _amplMax{};
291  double _eSum{};
292  double _energyLoss{};
293  double _clusterWidthX{},_clusterWidthY{};
294  double _ampl33{},_ampl55{},_ampl77{};
295  int _ncell33{},_ncell55{},_ncell77{};
296  int _storeNtuple{};
297  double _xLocalRecoCOG{},_yLocalRecoCOG{};
298  double _xLocalRecoEdge{},_yLocalRecoEdge{};
299  double _xLocalSim{},_yLocalSim{};
300 
301  std::vector <SimTrackerHitImplVec> _hitsInLayer{};
302 
303 
304 
305 
306 } ;
307 /** @} @} */
308 #endif
std::vector< TrackerHitImpl * > TrackerHitImplVec
Definition: CCDDigitizer.h:72
virtual Processor * newProcessor()
Definition: VTXDigitizer.h:136
std::vector< IonisationPoint > IonisationPointVec
Definition: CCDDigitizer.h:71
std::vector< SimTrackerHitImpl * > SimTrackerHitImplVec
Definition: CCDDigitizer.h:73
std::vector< SignalPoint > SignalPointVec
Definition: VTXDigitizer.h:41
std::vector< LCCollection * > LCCollectionVec
Definition: SiStripClus.h:55