PandaRoot
PndLmdDim.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 /*
14  * PndLmdDim.h
15  *
16  * this class gives you methods to retrieve parameters of the geometry
17  * dimensions are in cm and radian. this is a singleton and must be called
18  * with PndLmdDim::Instance()
19  *
20  * Global comments:
21  *
22  * TODO
23  *
24  *
25  * Created on: Oct 5, 2012
26  * Author: promme
27  *
28  * To use transformation functions you have to call Read_Transformation_matrices
29  *
30  * as of 2015 maintained by Roman Klasen, klasen@kph.uni-mainz.de
31 
32  */
33 
34 #ifndef PNDLMDDIM_H_
35 #define PNDLMDDIM_H_
36 
37 #include <cmath>
38 #include <TRandom.h>
39 #include <map>
40 #include <string>
41 #include <sstream>
42 #include <vector>
43 #include <iostream>
44 #include <TGeoManager.h>
45 #include <TGeoVolume.h>
46 #include <TGeoMatrix.h>
47 #include <TVector3.h>
48 #include <TMatrixT.h>
49 #include <TH2Poly.h>
50 #include <TMultiGraph.h>
51 #include <TGraph.h>
52 #include <stdlib.h>
53 #include <cctype>
54 
55 // work with DB
56 /* #include<PndLmdContFact.h> */
57 /* #include<TList.h> */
59 #include <TPolyLine.h>
60 /* #include "FairRuntimeDb.h" */
61 /* #include "FairRunAna.h" */
62 /* #include "FairRun.h" */
63 
64 typedef TMatrixT<double> TMatrixD; // hmm funny, should be declared in TMatrixT.h
65 
66 using namespace std;
67 
68 // simple class used for a fast particle track propagation
69 // based on a matrix fit to simulated data
70 // the representation of the basis is
71 // x, y, z, w, px/pz, py/pz, pw
72 // w und pw are weights in order to
73 // provide homogeneous coordinates
75  public:
76  double mom;
77  double M[7][7];
79  {
80  mom = 0;
81  for (int i = 0; i < 7; i++) {
82  for (int j = 0; j < 7; j++) {
83  M[i][j] = 0;
84  }
85  }
86  }
87  // set a matrix which is given as an array
88  void Set(const double *matrix)
89  {
90  for (int i = 0; i < 7 * 7; i++) {
91  M[i / 7][i % 7] = matrix[i];
92  }
93  }
94  // set a matrix with fitted hardcoded values
95  // up to now only 1.5 GeV/c and 15 GeV/c are available
96  // to produce more, please modify and use trafo_matrix_fit.C
97  PndLmdDimPropMat(double mom_init)
98  {
99  mom = mom_init;
100  if (mom == 1.5) {
101  double matrix_init[49] = {1.012, 1.22, 5.216e-05, 0.2599, 1.096, -0.1589, -0.0001459, -0.5851, 1.042, 0.001092,
102  -3.991e-06, 0.3103, 1.064, -3.991e-06, -0.04447, -0.04501, 0.001764, 11.24, -0.04385, 0.007135,
103  0.0001733, -1.36e-15, 1.816e-14, 3.031e-15, 1, 1.513e-16, -2.272e-16, 7.904e-17, 0.006668, 1.356,
104  -0.0001811, 0.000127, 0.9794, -0.1215, 0.4001, -0.6638, 0.04068, 0.001168, -7.627e-06, 0.2952,
105  0.9388, -7.627e-06, -1.003e-15, -3.248e-15, 9.487e-17, 5.594e-18, 2.41e-16, -1.086e-16, 1};
106  Set(matrix_init);
107  } else {
108  if (mom == 15) {
109  double matrix_init[49] = {0.9633, 0.2086, -0.01173, 0.2598, 1.123, -0.03217, -0.0002463, -0.104, 0.97, -0.006443,
110  5.119e-06, 0.06427, 1.12, 5.119e-06, -0.03828, 0.005853, 0.004478, 11.24, -0.0447, 0.001794,
111  2.613e-06, 2.376e-06, -0.0003806, -2.583e-07, 1, -1.647e-05, -1.871e-05, -3.862e-07, -0.01744, 0.1777,
112  -0.01384, 0.0001758, 0.9997, -0.02581, 0.4002, -0.09394, -0.05056, -0.002047, 3.624e-05, 0.06091,
113  0.994, 3.624e-05, 2.376e-06, -0.0003806, -2.583e-07, -3.862e-07, -1.647e-05, -1.871e-05, 1};
114  Set(matrix_init);
115  } else { // identity matrix
116  cout << " Warning in PndLmdDimPropMat: loading identity matrix only " << endl;
117  double matrix_init[49] = {1, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 1,
118  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 1};
119  Set(matrix_init);
120  }
121  }
122  }
123  // transform the position and momentum at the IP
124  // to a position and momentum at the first plane of the LMD
125  void Propagate(TVector3 &pos, TVector3 &momdir)
126  {
127  // TVector3 dir = mom.Unit();
128  double mommag = momdir.Mag();
129  // position transformation was evaluated for meter and angles * 10.
130  double xip[7] = {pos.X() / 100., pos.Y() / 100., pos.Z() / 100., 1, momdir.X() / momdir.Z() * 10., momdir.Y() / momdir.Z() * 10., 1};
131  double xlmd[7] = {0, 0, 0, 0, 0, 0, 0};
132  for (int i = 0; i < 7; i++) {
133  for (int j = 0; j < 7; j++) {
134  xlmd[i] += M[i][j] * xip[j];
135  }
136  }
137  pos.SetXYZ(xlmd[0] * 100., xlmd[1] * 100., xlmd[2] * 100.);
138  momdir.SetXYZ(xlmd[4] / 10., xlmd[5] / 10., 1.); // not an exact calculation ;)
139  momdir = momdir.Unit() * mommag;
140  }
141 };
142 
143 class Tkey {
144  public:
145  signed char half;
146  signed char plane;
147  signed char module;
148  signed char side;
149  signed char die;
150  signed char sensor;
151  bool operator<(const Tkey &comp) const
152  {
153  if (half < comp.half)
154  return true;
155  if (half > comp.half)
156  return false;
157  if (plane < comp.plane)
158  return true;
159  if (plane > comp.plane)
160  return false;
161  if (module < comp.module)
162  return true;
163  if (module > comp.module)
164  return false;
165  if (side < comp.side)
166  return true;
167  if (side > comp.side)
168  return false;
169  if (die < comp.die)
170  return true;
171  if (die > comp.die)
172  return false;
173  if (sensor < comp.sensor)
174  return true;
175  if (sensor >= comp.sensor)
176  return false;
177  return false;
178  }
179 
180  bool operator==(const Tkey &comp) const
181  {
182  return (half == comp.half) && (plane == comp.plane) && (module == comp.module) && (side == comp.side) && (die == comp.die) && (sensor == comp.sensor);
183  }
184 
185  Tkey(const Tkey &copy)
186  {
187  half = copy.half;
188  plane = copy.plane;
189  module = copy.module;
190  side = copy.side;
191  die = copy.die;
192  sensor = copy.sensor;
193  }
194 
195  Tkey(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
196  {
197  half = ihalf;
198  plane = iplane;
199  module = imodule;
200  side = iside;
201  die = idie;
202  sensor = isensor;
203  }
204 
205  Tkey(string key)
206  {
207  int sign(1);
208  int ikey(0); // 0,1,2,3,4,5 = ihalf, iplane, imodule, iside, idie, isensor
209  int number;
210  for (unsigned int ichar = 0; ichar < key.size(); ichar++) {
211  if (key[ichar] == '-') {
212  sign = -1;
213  }
214  if (isdigit(key[ichar])) {
215  number = (key[ichar] - '0') * sign;
216  sign = 1;
217  if (ikey == 0)
218  half = number;
219  if (ikey == 1)
220  plane = number;
221  if (ikey == 2)
222  module = number;
223  if (ikey == 3)
224  side = number;
225  if (ikey == 4)
226  die = number;
227  if (ikey == 5)
228  sensor = number;
229  ikey++;
230  }
231  }
232  if (ikey != 6)
233  cout << " Error in Generate_Tkey: key string " << key << " is not valid " << endl;
234  }
235  Tkey() {}
236 };
237 
238 class PndLmdDim {
239  private:
240  // in case you change anything in the geometry
241  // you must increase this version number in the corresponding .cxx file!
242  // it may affect the consistency between the geometry and the transformation matrices
243  static int geometry_version;
244 
245  static PndLmdDim *pinstance;
246  TGeoManager *fgGeoMan;
247  PndLmdDim();
248  PndLmdDim(const PndLmdDim &instance);
249  PndLmdDim &operator=(const PndLmdDim &instance)
250  {
251  this->box_size_x = instance.box_size_x; // to get rid from pedantic warnings
252  return *this;
253  }
254  ~PndLmdDim();
255  // the navigation paths for the detector geometry are stored here
256  vector<string> nav_paths;
257 
258  public:
259  static PndLmdDim &Get_instance();
260  static PndLmdDim *Instance();
261 
262  // this offset in the sensor id is introduced by
263  // pandaroot when assembling several detector components
264  // the lmd alone would have a sensor id range between 0 and n total sensors
265  // A detector loaded before the lmd will introduce an offset in
266  // that ID which can be determined from a loaded root geometry by calling
267  // Set_sensIDoffset(-1);
268  // TODO: add everywhere that sensID offset
269  unsigned int sensIDoffset;
270 
271  // set the sensIDoffset
272  // in case offset is < 0 the offset is
273  // tried to be determined from a possibly loaded root geometry
274  bool Set_sensIDoffset(int offset = -1);
275 
276  // FairRun* ana;
277  // FairRuntimeDb* rtdb;
278  // pi
279  double pi;
280  // number of detector planes
281  unsigned int n_planes;
282  // number of modules per plane half
283  unsigned int nmodules;
284  // position of planes where the first plane defines the origin
285  double *plane_pos_z;
286  // ****************************** plane half array *****************************
287  // x in lmd reference frame
289  // y in lmd reference frame
291  // z in lmd reference frame
293  // phi rotation in lmd reference frame
295  // theta rotation in lmd reference frame
297  // psi rotation in lmd reference frame
299  // ****************************** plane half supports ******************************
300  // x in the detector half reference frame
302  // y in the detector half reference frame
304  // z in the detector half reference frame
306  // phi rotation in the detector half reference frame
308  // theta rotation in the detector half reference frame
310  // psi rotation in the detector half reference frame
312 
313  // ****************************** cvd cooling support discs ************************
314  // cvd_diamond is cut out of 79.5 mm discs of 200 micron thickness
315  // inner min. radius due to beam pipe + a safety margin
316  double inner_rad;
317  // not used yet but should be the outer acceptance
318  double outer_rad;
319  // number of CVD diamond discs per plane
321  // radius of a CVD diamond disc
322  double cvd_disc_rad;
323  // the half of the diamond thickness
325  // even and odd discs in a plane will be shifted in z in order to prevent
326  // mechanical damage during assembly
328  // angle from the division of a circle into n_cvd_discs
329  double delta_phi;
330  // a polygon of n_cvd_discs sides fitting a radius of inner_rad
331  // has a side length pol_side_lg of
333  // the minimum distance to the center of the polygone is given by
335  // the cvd disc has to be placed such that the disc crosses
336  // the inner ring at an angle of 0 and delta_phi
337  // this defines the distance to the center according to pythagoras
339 
340  // the mechanical alignment precision is defined as an offset of and tilt around
341  // the middle of the cvd diamond.
342  // Values are standard deviation.
343  // first comes translation than rotation
344 
345  // x is radial to the beam pipe when support is aligned
346  double cvd_offset_x;
347  // y is tangent to the beam pipe when support is aligned
348  double cvd_offset_y;
349  // z is along the beam pipe when support is aligned
350  double cvd_offset_z;
351  // phi rotation in the reference frame of the plane half support
352  double cvd_tilt_phi;
353  // theta rotation in the reference frame of the plane half support
355  // psi rotation in the reference frame of the plane half support
356  double cvd_tilt_psi;
357  // thickness of the kapton flexible circuits to the sensors
358  // (chosen to be thicker to simulate also the influence of the
359  // printed circuits themselves)
361  // *********************************** one side on a CVD disc *************************************
362  // x is radial to the beam pipe when cvd is aligned
364  // y is tangent to the beam pipe when cvd is aligned
366  // z is along the beam pipe when cvd is aligned
367  double side_offset_z; // should not be used due to clashing volumes with the support!
368  // phi rotation in the reference frame of the cvd support
369  double side_tilt_phi; // should not be used due to clashing volumes with the support!
370  // theta rotation in the reference frame of the cvd support
372  // psi rotation in the reference frame of the cvd support
373  double side_tilt_psi; // should not be used due to clashing volumes with the support!
374  // *********************************** HV-MAPS *************************************
375  //
376  // left right
377  //
378  // |---------|----------|----------| top
379  // ||------|-||-------|-||-------|-|
380  // || | || | || | |
381  // || 1 | || 2 | || 3 | | row 1
382  // || | || active| || | |
383  // ||------|-||-------|-||-------|-|
384  // | | passive | |
385  // |---------|----------|----------| bottom
386  // gap
387  // |----------|----------| top
388  // | passive | |
389  // ||-------|-||-------|-|
390  // || | || | |
391  // || 2 | || 3 | | row 2
392  // || active| || | |
393  // ||-------|-||-------|-|
394  // |----------|----------| bottom
395  //
396  // left right
397  //
398  // A y; maps_n_row; height
399  // |
400  // |
401  // --> x; maps_n_col; width
402  //
403  // the current design foresees a rotation of the first parameters
404 
405  // even when several maps are placed on one die
406  // those will be placed in the simulation next to
407  // each other as separate detectors
408  // assumption: sensor 1 and 2 are on one die, the rest is separated
411  // enabled [row][col]
412  bool **enabled;
413  // number of sensors per side
415  // NOTE: MOST of the following VARIABLES are HALF of it
416  // due to geometry construction in GEANT
424 
425  double maps_active_pixel_size; // pixel size NOT half of it
426 
427  double maps_width;
428  double maps_height;
429 
432 
433  double die_gap; // (cm)
434 
437 
438  // the mechanical alignment precision is defined as an offset of and tilt around
439  // the middle of the cvd diamond.
440  // Values are standard deviation.
441  // first comes translation than rotation
442  // translations along z as well as rotations around x and y are
443  // negligible for dies glued on a cvd diamond
444  // rotation around z is not working yet
445 
446  // x is along the edge of the cvd disc
447  double die_offset_x;
448  // y is orthogonal to the edge of the cvd disc
449  double die_offset_y;
450  // z is along the beam pipe
451  double die_offset_z; // should not be used -> crashing volumes;
452  // x is a rotation around the edge of the cvd disc
453  double die_tilt_phi; // please do not use yet
454  // y is a rotation around the orthogonal component of the edge of the cvd disc
455  double die_tilt_theta; // please do not use yet
456  // z is a rotation around an axis parallel to the along the beam pipe
457  double die_tilt_psi; // please do not use yet;
458  // *********************************** one sensor *************************************
459  // x is mostly radial to the beam pipe
461  // y is mostly tangent to the beam pipe
463  // z is along the beam pipe
464  double sensor_offset_z; // should not be used due to clashing volumes with cvd and flex prints
465  // phi rotation in the reference frame of the sensor
467  // theta rotation in the reference frame of the sensor
468  double sensor_tilt_theta; // should not be used due to clashing volumes with cvd and flex prints
469  // psi rotation in the reference frame of the sensor
470  double sensor_tilt_psi; // should not be used due to clashing volumes with cvd and flex prints
471  //*********************************** lumi box parameters ***********************************
472  // see CAD files for details
473  // https://edms.cern.ch/nav/P:FAIR-000000719:V0/P:FAIR-000000726:V0/TAB3
474  // width
475  double box_size_x;
476  // height
477  double box_size_y;
478  // length
479  double box_size_z;
480  // thickness of the V2A steel plates
482  // position of the inner rib
483  double pos_rib;
484  // beam pipe radius at entrance
485  double rad_entrance;
486  // beam pipe radius at exit
487  double rad_exit;
488  // beam pipe separating non interacting paricles
489  double rad_pipe;
490  // beam pipe thickness;
492  // cone height of the transition region
494  // length of the inner pipe
495  double length_pipe;
496  // position of the first detector plane
497  double pos_plane_0;
498  //*********************************** global parameters *************************************
499  // where bending starts with
501  // the bending radius
502  double r_bend;
503  // and the angle of the circle path
504  double phi_bend;
505  // the point where both tangents of the straight beam pipe tubes meet is
506  double pos_rot_z;
507  // z position of the lmd
508  double pos_z; //(cm)
509  double end_seg_bend;
510  // x position of the lmd
511  double pos_x;
512  // y position of the lmd
513  double pos_y;
514  double rot_x;
515  double rot_y;
516  double rot_z;
517 
518  // returns false when one of the variables exceeds design values
519  bool Is_valid_idcall(int ihalf, int iplane = 0, int imodule = 0, int iside = 0, int idie = 0, int isensor = 0)
520  {
521  if (ihalf < 0 || ihalf >= 2)
522  return false;
523  if (iplane < 0 || (unsigned)iplane >= n_planes)
524  return false;
525  if (imodule < 0 || imodule >= n_cvd_discs)
526  return false;
527  if (iside < 0 || iside >= 2)
528  return false;
529  if (idie < 0 || idie >= 2)
530  return false;
531  // allow to count the non existing inner sensor at die 2
532  if (isensor < 0 || isensor >= n_sensors + 1)
533  return false;
534  return true;
535  }
536 
537  // get the sensor id for a sensor on a given side, module and plane
538  int Get_sensor_id(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
539  {
540  if (idie == 1)
541  isensor += 2; // the parallel sensor to sensor 0 is not there!
542  int result = isensor + (iside + (imodule + (iplane + ihalf * n_planes) * nmodules) * 2) * n_sensors;
543  return result + sensIDoffset;
544  }
545 
546  // get the sensor position by it's id in terms of plane module and side
547  void Get_sensor_by_id(const int sensor_id, int &ihalf, int &iplane, int &imodule, int &iside, int &idie, int &isensor)
548  {
549  int _sensor_id = sensor_id - sensIDoffset;
550  isensor = _sensor_id % n_sensors;
551  idie = 0;
552  if (isensor > 2) {
553  idie = 1;
554  isensor -= 2; // add the first but non existing sensor at die 2
555  }
556  _sensor_id /= n_sensors;
557  iside = _sensor_id % 2;
558  _sensor_id /= 2;
559  imodule = _sensor_id % nmodules;
560  _sensor_id /= nmodules;
561  iplane = _sensor_id % n_planes;
562  _sensor_id /= n_planes;
563  ihalf = _sensor_id % 2;
564  if (!Is_valid_idcall(ihalf, iplane, imodule, iside, idie, isensor)) {
565  ihalf = 0;
566  iplane = 0;
567  imodule = 0;
568  iside = 0;
569  idie = 0;
570  isensor = 0;
571  cout << "Error in PndLmdDim::Get_sensor_by_id: " << sensor_id << " is not a valid sensor id!" << endl;
572  }
573  }
574 
575  // when calling for a displaced module, the generated
576  // value must be kept same for the cvd diamond as well
577  // as the sensors sitting on the diamond it self
578  // the storage is realized by a map
579  // the vector contains offsets in x, y, z, rotphi, rottheta, rotpsi;
580  map<Tkey, vector<double>> offsets;
581  map<Tkey, vector<double>>::iterator itoffset;
582  /* replaced by an ugly but faster version below
583  string Generate_key(int ihalf, int iplane, int imodule, int iside, int idie, int isensor){
584  stringstream keystream;
585  keystream << ihalf << iplane << imodule << iside << idie << isensor;
586  return keystream.str();
587  }*/
588 
594  char *itoa(int value, char *result, int base)
595  {
596  // check that the base if valid
597  char *last_char;
598  if (base < 2 || base > 36) {
599  *result = '\0';
600  return result;
601  }
602  char *ptr = result, *ptr1 = result, tmp_char;
603  int tmp_value;
604  do {
605  tmp_value = value;
606  value /= base;
607  *ptr++ = "zyxwvutsrqponmlkjihgfedcba9876543210123456789abcdefghijklmnopqrstuvwxyz"[35 + (tmp_value - value * base)];
608  } while (value);
609  // Apply negative sign
610  if (tmp_value < 0)
611  *ptr++ = '-';
612  last_char = ptr;
613  *ptr-- = '\0';
614  // cout << last_char << endl;
615  while (ptr1 < ptr) {
616  tmp_char = *ptr;
617  *ptr-- = *ptr1;
618  *ptr1++ = tmp_char;
619  }
620  return last_char;
621  }
622 
623  string Generate_key(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
624  {
625  char key[100];
626  char *ptr;
627  ptr = itoa(ihalf, key, 10);
628  ptr = itoa(iplane, ptr, 10);
629  ptr = itoa(imodule, ptr, 10);
630  ptr = itoa(iside, ptr, 10);
631  ptr = itoa(idie, ptr, 10);
632  ptr = itoa(isensor, ptr, 10);
633  string result(key);
634  // stringstream keystream;
635  // keystream << ihalf << iplane << imodule << iside << idie << isensor;
636  return result;
637  }
638 
639  // generate a unique integer key not same as the string above since there negative numbers are allowed
640  // so in case of adding -1's key is not unique!
641  // moreover numbers should be kept 1 digit long
642  int Generate_keynumber(unsigned int ihalf = 0, unsigned int iplane = 0, unsigned int imodule = 0, unsigned int iside = 0, unsigned int idie = 0, unsigned int isensor = 0)
643  {
644  stringstream keystream;
645  keystream << ihalf << iplane << imodule << iside << idie << isensor;
646  int key;
647  key = atoi(keystream.str().c_str());
648  return key;
649  }
650 
651  // same structure as for offsets is used for the transformation matrices
652  // stored are matrix operations
653  // global -> local lumi
654  // key: ihalf = -1 iplane = -1 imodule = -1 iside = -1 idie = -1 isensor = -1
655  // local lumi -> local side on cvd disc
656  // key: ihalf >=0 iplane >= 0 imodule >=0 iside = -1 idie = -1 isensor = -1
657  // local side on cvd disc -> local sensor
658  // key: all variable
659  // map<string, TGeoMatrix* > transformation_matrices;
660  // map<string, TGeoMatrix* > transformation_matrices_aligned; // alternative aligned detector description
661  // map<string, TGeoMatrix* >::iterator it_transformation_matrices;
662  // to increase the performance by a factor of 5-6 a struct is used for the key now
663  map<Tkey, TGeoMatrix *> transformation_matrices;
664  map<Tkey, TGeoMatrix *> transformation_matrices_aligned; // alternative aligned detector description
665  map<Tkey, TGeoMatrix *>::iterator it_transformation_matrices;
666 
667  // cleanup some maps containing only references
668  void Cleanup();
669 
670  // read transformation matrices from a given file
671  // aligned and not aligned are two separate maps
672  // containing the description of the detector positions
673  // if not filename is specified matrices are searched in
674  // VMCWORKDIR/input/matrices.txt
675  // you may overwrite the version number if necessary
676  // VMCWORKDIR/geometry folder is used if no filename is specified
677  // ATTENTION! aligned matrices means perfect geometry, set to false if using misaligned geometry
678  void Read_transformation_matrices(string filename = "", bool aligned = true, int version_number = geometry_version);
679 
680  // write transformation matrices from a given file
681  // aligned and not aligned are two separate maps
682  // containing the description of the detector positions
683  // you may overwrite the version number if necessary
684  // version == 0 is reserved for backward compatibility!
685  // VMCWORKDIR/geometry folder is used if no filename (e.g. "") is specified
686  // warning overwrites existing trafo_matrices_lmd.dat!
687  void Write_transformation_matrices(string filename, bool aligned = true, int version_number = geometry_version);
688 
689  // read transformation matrices from a loaded geometry
690  // aligned and not aligned are two separate maps
691  // containing the description of the detector positions
692  // the geometry must be loaded otherwise matrices cannot be read
693  // version number will be set according to the geometry version number
694  // ToDo: multiply also matrices on the way to the key matrices
695  // in case those are not unity matrices
696  bool Read_transformation_matrices_from_geometry(bool aligned = true);
697 
698  // apply transformation matrices to a loaded geometry
699  // aligned and not aligned are two separate maps
700  // containing the description of the detector positions
701  // the geometry must be loaded otherwise matrices cannot be read
702  // version number will be set according to the geometry version number
703  // IMPORTANT: you may choose which PndLmdDim matrices you want to use
704  // but a ROOT Geometry can always be only aligned. The original
705  // matrix stays untouched!
706  // TODO: multiply also matrices on the way to the key matrices
707  // in case those are not unity matrices
708  // TODO: Find out how to store the aligned geometry as a default
709  // one to pandaroot parameter files
710  bool Write_transformation_matrices_to_geometry(bool aligned = true);
711 
712  // get a list of sensor paths in the geometry navigation model
713  // returns the path to the lmd top volume
714  // It is a recursive search, call it once with the default found_lmd variable
715  // in case first_call then gGeoManager->CdTop(); is executed to get to the top node
716  // of a geometry
717  // The geometry must be loaded
718  string Get_List_of_Sensors(vector<string> &list_of_sensors, bool found_lmd = false, bool first_call = true);
719 
720  // check a list of sensor paths for validity
721  // result is true if tests were sucessful
722  // offset is the offset in the sensor number which may be
723  // not 0 if the geometry was not created in the first place
724  bool Test_List_of_Sensors(vector<string> list_of_sensors, int &offset);
725 
726  // Get an offset for a volume, if not existent and random
727  // a random offset is generated and stored
728  // Is used during generation of geometries when calling
729  void Get_offset(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, double &x, double &y, double &z, double &rotphi, double &rottheta, double &rotpsi,
730  bool random = false);
731 
732  // set an offset for example for the case of existent alignment
733  // values, to generate with Get_offset a geometry that matches
734  // those offsets
735  void Set_offset(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, double x, double y, double z, double rotphi, double rottheta, double rotpsi);
736 
737  // read alignment constants from DB
738  // Feb 2013: currently DB is ASCII file
739  void Read_DB_offsets(PndLmdAlignPar *lmdalignpar);
740 
741  // correct transformation matrices by align
742  void Correct_transformation_matrices();
743  void reCreate_transformation_matrices();
744  // several functions returning the position and orientation of
745  // the luminosity detector
746 
747  // the global system is defined by PANDA
748  // returns the position and rotation of the lmd detector
749  // the origin is PANDA interaction point
750  // first comes translation than rotation
751  void Get_pos_lmd_global(double &x, double &y, double &z, double &rotx, double &roty, double &rotz, bool misaligned = false)
752  {
753  rotx = rot_x;
754  roty = rot_y; // phi_bend;
755  rotz = rot_z;
756  x = pos_x;
757  y = pos_y;
758  z = pos_z;
759  if (misaligned)
760  rotx += 0; // pedantic compiler fix
761  }
762 
763  // decoding a digital hit by the position of the sensor given by the sensor ID
764  // and the row and column entry
765  // parameters of MC geometry input are used to determine the position of
766  // the active area within the sensor volume
767  // The hit point is returned in the panda global reference frame
768  // needs transformation matrices to determine the position
769  // column and row can be also the mean from a cluster and therefore
770  // not an integer
771  TVector3 Decode_hit(const int sensorID, const double column, const double row, const bool aligned = true, bool newVersion = false);
772 
773  // matrices for a fast prediction of a particle track to the lmd
774  // depending on the momentum setting of the hesr which is the key
775  map<double, PndLmdDimPropMat> propagation_matrices;
776 
777  // propagate a particle at the IP to the first plane of the lmd
778  // the propagation is based on a transformation matrix which was
779  // fit to GEANT4 propagated data from January 2015
780  // only 1.5 GeV/c and 15 GeV/c are implemented properly
781  // TODO: interpolation between matrices
782  // back propagation not implemented yet since
783  // Mathematica inverted matrices did not work
784  void Propagate_fast_ip_to_lmd(TVector3 &pos, TVector3 &mom, double pbeam);
785 
786  // Find the corresponding sensor on the opposite side of a module
787  // point is in the panda frame
788  // point is transformed into the local lmd frame and
789  // the projection is tested 1. to lie on the given sensor
790  // 2. on all other sensors on the opposite side for intersection
791  // returns true if the corresponding sensor was found
792  // in that case input parameters are set to the corresponding sensor
793  // otherwise those are set to -1!
794  // do not use in time critical cases!
795  bool Get_overlapping_sensor(const TVector3 &point, int &ihalf, int &iplane, int &imodule, int &iside, int &idie, int &isensor, bool aligned = true);
796 
797  // Test if a point lies on a specific sensor
798  // The point in the panda frame is transformed into the sensor frame
799  // the projection is tested to be within the boundaries of the active sensor area
800  bool Is_on_Sensor(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
801 
802  // get the overlapping sensors by geometric constraints
803  // hard coded overlapping areas are given sorted by the overlapping area
804  // returned is the number of overlapping sensors
805  // jdie and jsensor are vectors with the length of the returned number
806  int Get_overlapping_sensor(int idie, int isensor, vector<int> &jdie, vector<int> &jsensor);
807 
808  // get the transformation matrix from the PANDA global reference frame to the
809  // Luminosity reference frame
810  TGeoHMatrix Get_transformation_global_to_lmd_local(bool aligned = true);
811 
812  // get the transformation matrix from lmd local reference frame to the
813  // cvd disc surface reference frame
814  TGeoHMatrix Get_transformation_lmd_local_to_module_side(int ihalf, int iplane, int imodule, int iside, bool aligned = true);
815 
816  // get the transformation matrix from lmd cvd disc surface frame to the
817  // sensor reference frame
818  TGeoHMatrix Get_transformation_module_side_to_sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
819 
820  // get the transformation matrix from PANDA global reference frame to the
821  // sensor reference frame
822  TGeoHMatrix Get_transformation_global_to_sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
823 
824  // get the transformation matrix from lumi reference frame to the
825  // sensor reference frame
826  TGeoHMatrix Get_transformation_lmd_local_to_sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
827 
828  // get the inverse transformation matrix from the PANDA global reference frame to the
829  // Luminosity reference frame
830  TGeoHMatrix Get_transformation_lmd_local_to_global(bool aligned = true);
831 
832  // get the inverse transformation matrix from lmd local reference frame to the
833  // cvd disc surface reference frame
834  TGeoHMatrix Get_transformation_module_side_to_lmd_local(int ihalf, int iplane, int imodule, int iside, bool aligned = true);
835 
836  // get the inverse transformation matrix from lmd cvd disc surface frame to the
837  // sensor reference frame
838  TGeoHMatrix Get_transformation_sensor_to_module_side(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
839 
840  // get the inverse transformation matrix from PANDA global reference frame to the
841  // sensor reference frame
842  TGeoHMatrix Get_transformation_sensor_to_global(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
843 
844  // get the inverse transformation matrix from lmd local reference frame to the
845  // sensor reference frame
846  TGeoHMatrix Get_transformation_sensor_to_lmd_local(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
847 
848  // get the transformation matrix from a ideal sensor to the aligned one
849  TGeoHMatrix Get_transformation_sensor_to_sensor_aligned(int ihalf, int iplane, int imodule, int iside, int idie, int isensor);
850 
851  // get the transformation matrix from a aligned sensor to the ideal one
852  TGeoHMatrix Get_transformation_sensor_aligned_to_sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor);
853 
854  //************************** 3d vector transformations ******************************
855 
856  // Transform from the PANDA global reference frame to the
857  // Luminosity reference frame
858  TVector3 Transform_global_to_lmd_local(const TVector3 &point, bool isvector = false, bool aligned = true);
859 
860  // Transform from lmd local reference frame to the
861  // cvd disc surface reference frame
862  TVector3 Transform_lmd_local_to_module_side(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, bool isvector = false, bool aligned = true);
863 
864  // Transform from lmd cvd disc surface frame to the
865  // sensor reference frame
866  TVector3 Transform_module_side_to_sensor(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
867 
868  // Transform from PANDA global reference frame to the
869  // sensor reference frame
870  TVector3 Transform_global_to_sensor(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
871 
872  // Transform from lmd local reference frame to the
873  // sensor reference frame
874  TVector3 Transform_lmd_local_to_sensor(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
875 
876  // Transform from the PANDA global reference frame to the
877  // Luminosity reference frame
878  TVector3 Transform_lmd_local_to_global(const TVector3 &point, bool isvector = false, bool aligned = true);
879 
880  // Transform from lmd local reference frame to the
881  // cvd disc surface reference frame
882  TVector3 Transform_module_side_to_lmd_local(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, bool isvector = false, bool aligned = true);
883 
884  // Transform from lmd cvd disc surface frame to the
885  // sensor reference frame
886  TVector3 Transform_sensor_to_module_side(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
887 
888  // Transform from PANDA global reference frame to the
889  // sensor reference frame
890  TVector3 Transform_sensor_to_global(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
891 
892  // Transform from lmd local reference frame to the
893  // sensor reference frame
894  TVector3 Transform_sensor_to_lmd_local(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false, bool aligned = true);
895 
896  // Transform from a ideal sensor to the aligned one
897  TVector3 Transform_sensor_to_sensor_aligned(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false);
898 
899  // Transform from a aligned sensor to the ideal one
900  TVector3 Transform_sensor_aligned_to_sensor(const TVector3 &point, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool isvector = false);
901 
902  // Transform from one sensor to an other
903  TVector3 Transform_sensor_to_sensor(const TVector3 &point, int ihalf_from, int iplane_from, int imodule_from, int iside_from, int idie_from, int isensor_from, int ihalf_to,
904  int iplane_to, int imodule_to, int iside_to, int idie_to, int isensor_to, bool isvector = false, bool aligned = true);
905 
906  // ************************* 3d matrix rotations ***************************************
907 
908  // Transform from the PANDA global reference frame to the
909  // Luminosity reference frame
910  // treats only 3 x 3 matrices representing the space
911  TMatrixD Transform_global_to_lmd_local(const TMatrixD &matrix, bool aligned = true);
912 
913  // Transform from lmd local reference frame to the
914  // cvd disc surface reference frame
915  // treats only 3 x 3 matrices representing the space
916  TMatrixD Transform_lmd_local_to_module_side(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, bool aligned = true);
917 
918  // Transform from lmd local reference frame to the
919  // sensor reference frame
920  // treats only 3 x 3 matrices representing the space
921  TMatrixD Transform_lmd_local_to_sensor(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
922 
923  // Transform from lmd cvd disc surface frame to the
924  // sensor reference frame
925  // treats only 3 x 3 matrices representing the space
926  TMatrixD Transform_module_side_to_sensor(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
927 
928  // Transform from PANDA global reference frame to the
929  // sensor reference frame
930  // treats only 3 x 3 matrices representing the space
931  TMatrixD Transform_global_to_sensor(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
932 
933  // Transform from the PANDA global reference frame to the
934  // Luminosity reference frame
935  // treats only 3 x 3 matrices representing the space
936  TMatrixD Transform_lmd_local_to_global(const TMatrixD &matrix, bool aligned = true);
937 
938  // Transform from lmd local reference frame to the
939  // cvd disc surface reference frame
940  // treats only 3 x 3 matrices representing the space
941  TMatrixD Transform_module_side_to_lmd_local(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, bool aligned = true);
942 
943  // Transform from lmd cvd disc surface frame to the
944  // sensor reference frame
945  // treats only 3 x 3 matrices representing the space
946  TMatrixD Transform_sensor_to_module_side(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
947 
948  // Transform from sensor reference frame to the
949  // lmd local reference frame
950  // treats only 3 x 3 matrices representing the space
951  TMatrixD Transform_sensor_to_lmd_local(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
952 
953  // Transform from PANDA global reference frame to the
954  // sensor reference frame
955  // treats only 3 x 3 matrices representing the space
956  TMatrixD Transform_sensor_to_global(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
957 
958  // Transform from a ideal sensor to the aligned one
959  // treats only 3 x 3 matrices representing the space
960  TMatrixD Transform_sensor_to_sensor_aligned(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor);
961 
962  // Transform from a aligned sensor to the ideal one
963  // treats only 3 x 3 matrices representing the space
964  TMatrixD Transform_sensor_aligned_to_sensor(const TMatrixD &matrix, int ihalf, int iplane, int imodule, int iside, int idie, int isensor);
965 
966  // Transform from one sensor to an other
967  // treats only 3 x 3 matrices representing the space
968  TMatrixD Transform_sensor_to_sensor(const TMatrixD &matrix, int ihalf_from, int iplane_from, int imodule_from, int iside_from, int idie_from, int isensor_from, int ihalf_to,
969  int iplane_to, int imodule_to, int iside_to, int idie_to, int isensor_to, bool aligned = true);
970 
971  // *************************
972 
973  // get a pointer to the requested matrices with checks
974  // returns nullptr if no matrices available
975  // do not delete!
976  map<Tkey, TGeoMatrix *> *Get_matrices(bool aligned = true);
977 
978  // get a pointer to the requested matrix with checks
979  // returns nullptr if no matrix available
980  // do not delete or modify unless you know why you want to
981  // load the two files for matrices first
982  TGeoMatrix *Get_matrix(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
983 
984  // get the transformation matrix for an path in an existing root geometry
985  // no checks are performed in advance
986  // result may be 0!
987  // if (aligned) the matrix after a possible alignment is returned
988  // in that case details to the matrix must be provided in form
989  // of ihalf ... isensor
990  // TODO: get rid of path and do it only on the basis of ihalf ... isensor
991  // if (!aligned) the original matrix is returned
992  TGeoHMatrix *Get_matrix(string path, bool aligned = true, int ihalf = -1, int iplane = -1, int imodule = -1, int iside = -1, int idie = -1, int isensor = -1);
993 
994  // set the transformation matrix for an path in an existing root geometry
995  // A matrix can be only aligned,
996  // therefore by default original matrices are not touched!
997  // since a key must be created for a node
998  // details to it must be provided in form of
999  // ihalf ... isensor
1000  // TODO: get rid of path and do it only on the basis of ihalf ... isensor
1001  bool Set_matrix(string path, TGeoHMatrix *matrix, int ihalf = -1, int iplane = -1, int imodule = -1, int iside = -1, int idie = -1, int isensor = -1); //, bool aligned = true);
1002 
1003  // get the difference between two matrices of the aligned and misaligned branch
1004  // in terms of displacement and euler angles
1005  // all 0 in case of troubles and result is false
1006  // load the two files for matrices first
1007  // NOTE: output values are formatted to µm and µrad!!!
1008  bool Get_matrix_difference(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, double &dx, double &dy, double &dz, double &dphi, double &dtheta, double &dpsi);
1009 
1010  // calculates the differences of matrices between loaded
1011  // aligned and misaligned ones
1012  // As an output a table and some histograms are generated
1013  // load the two files for matrices first
1014  void Calc_matrix_offsets();
1015 
1016  // see get matrix
1017  TGeoMatrix *Get_matrix_global_to_lmd_local(bool aligned = true);
1018 
1019  // see get matrix
1020  TGeoMatrix *Get_matrix_lmd_local_to_module_side(int ihalf, int iplane, int imodule, int iside, bool aligned = true);
1021 
1022  // see get matrix
1023  TGeoMatrix *Get_matrix_module_side_to_sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true);
1024 
1025  // TGeoMatrix* Get_matrix_global_to_sector_local(int ihalf, int iplane, int imodule, bool aligned=true);
1026 
1027  // x, y, z coordinate transformation from the PANDA global reference frame to the
1028  // local reference frame of the luminosity monitor
1029  void Transform_global_to_lmd_local(double &x, double &y, double &z, bool aligned = true);
1030 
1031  // x, y, z vector transformation from the PANDA global reference frame to the
1032  // local reference frame of the luminosity monitor
1033  void Transform_global_to_lmd_local_vect(double &x, double &y, double &z, bool aligned = true);
1034 
1035  // x, y, z coordinate transformation from the local sensor reference frame to the
1036  // local reference frame of the luminosity monitor
1037  // void transform_sensor_local_to_lmd_local(const int sensor_id, double& x, double& y, double& z, bool misaligned = false);
1038 
1039  // x, y, z coordinate transformation from the local luminosity reference frame to the
1040  // global reference frame of panda
1041  // void transform_local_lmd_to_global(const int sensor_id, double& x, double& y, double& z, bool misaligned = false);
1042 
1043  // x, y, z coordinate transformation between the local sensor reference frames
1044  // of misaligned and aligned sensors
1045  // void transform_local_sensor();
1046 
1047  // Generates the luminosity monitor geometry into the mother volume
1048  // Please make sure that mother volume is large enough or that it
1049  // is an assembly volume
1050  // in addition rotation and translation matrices are calculated
1051  // to store those into a file please use
1052  // Write_transformation_matrices(filename, false);
1053  // in case of misaligned, offsets are retrieved via get_offset
1054  // and included into the geometry
1055  void Generate_rootgeom(TGeoVolume &mothervol, bool misaligned = false);
1056 
1057  // returns false
1058  // if version number of the geometry could not be retrieved from a loaded geometry
1059  // geometry must be available via the root geometry manager
1060  bool Retrieve_version_number();
1061 
1062  // small function to test some transformation matrices and methods
1063  void Test_matrices();
1064 
1065  // Draw the Sensors as overlays to an active root pad
1066  // Projection in XY is used
1067  // lmd_frame == true : The Lumi reference frame is used
1068  void Draw_Sensors(int iplane, bool aligned = true, bool lmd_frame = true, int glside = 2); // by default draw both sides
1069  // Get one Sensor as a polyline to be drawn to an active root pad
1070  // Projection in XY is used
1071  // Important -> Do not delete the PolyLine until you are sure that you are
1072  // finished with displaying it!
1073  // lmd_frame == true : The Lumi reference frame is used
1074  TPolyLine *Get_Sensor_Shape(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true, bool lmd_frame = true);
1075 
1076  // Get one Sensor as a vector of graphs
1077  // Projection in XY is used
1078  // lmd_frame == true : The Lumi reference frame is used
1079  // else the sensor reference frame
1080  // if pixel_subdivision: several graphs are returned
1081  // representing the passive area and the
1082  // active area subdivided into several pixels
1083  // please do not call several times with same
1084  // parameters if previous result and it's contents
1085  // are not deleted
1086  vector<TGraph *>
1087  Get_Sensor_Graph(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true, bool lmd_frame = true, bool pixel_subdivision = true);
1088 
1089  // Get a TH2Poly histogram where bins matching the x_y projected shapes
1090  // of the sensors in one plane
1091  // The active area can be subdivided further
1092  // the active pixels: NOT RECOMMENDED due to heavy ram usage
1093  // Please rename it if you intend to call the function
1094  // several times with same parameters
1095  TH2Poly *Get_histogram_Plane(int iplane, int iside, bool aligned = true, bool lmd_frame = true, bool pixel_subdivision = false);
1096 
1097  // Get a TH2Poly histogram for one module side
1098  // in the reference frame of the lumi
1099  // on request in the panda reference frame
1100  // The active area can be subdivided further
1101  // into the active pixels: WARNING heavy ram usage
1102  // Please rename it if you intend to call the function
1103  // several times with same parameters
1104  TH2Poly *Get_histogram_Moduleside(int ihalf, int iplane, int imodule, int iside, bool aligned = true, bool lmd_frame = true, bool pixel_subdivision = true);
1105 
1106  // Get a TH2Poly histogram for one sensor
1107  // in the reference frame of the lumi
1108  // on request in the panda reference frame
1109  // The active area can be subdivided into
1110  // and the active pixels
1111  // Please rename it if you intend to call the function
1112  // several times with same parameters
1113  TH2Poly *Get_histogram_Sensor(int ihalf, int iplane, int imodule, int iside, int idie, int isensor, bool aligned = true, bool lmd_frame = true);
1114 
1115  // gets all ids that correspond to an overlapping area
1116  std::vector<int> getAvailableOverlapIDs();
1117  int makeOverlapID(int firstSensorId, int secondSensorId);
1118  int getID1fromOverlapID(int overlapID);
1119  int getID2fromOverlapID(int overlapID);
1120 
1121  // this is a legacy function and should not be called anymore
1122  int makeModuleID(int overlapID);
1123 };
1124 
1125 #endif /* PNDLMDDIM_H_ */
double maps_active_pixel_size
Definition: PndLmdDim.h:425
double sensor_tilt_psi
Definition: PndLmdDim.h:470
double rot_z
Definition: PndLmdDim.h:516
double die_tilt_psi
Definition: PndLmdDim.h:457
double pos_y
Definition: PndLmdDim.h:513
double r_bend
Definition: PndLmdDim.h:502
map< Tkey, vector< double > >::iterator itoffset
Definition: PndLmdDim.h:581
void Get_pos_lmd_global(double &x, double &y, double &z, double &rotx, double &roty, double &rotz, bool misaligned=false)
Definition: PndLmdDim.h:751
Tkey(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
Definition: PndLmdDim.h:195
double half_tilt_phi
Definition: PndLmdDim.h:294
char * itoa(int value, char *result, int base)
Definition: PndLmdDim.h:594
bool Is_valid_idcall(int ihalf, int iplane=0, int imodule=0, int iside=0, int idie=0, int isensor=0)
Definition: PndLmdDim.h:519
double plane_half_tilt_theta
Definition: PndLmdDim.h:309
string Generate_key(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
Definition: PndLmdDim.h:623
double half_offset_z
Definition: PndLmdDim.h:292
double pos_z
Definition: PndLmdDim.h:508
double plane_half_offset_y
Definition: PndLmdDim.h:303
map< Tkey, TGeoMatrix * >::iterator it_transformation_matrices
Definition: PndLmdDim.h:665
double cvd_offset_y
Definition: PndLmdDim.h:348
double die_tilt_theta
Definition: PndLmdDim.h:455
double die_offset_y
Definition: PndLmdDim.h:449
double delta_phi
Definition: PndLmdDim.h:329
double die_tilt_phi
Definition: PndLmdDim.h:453
map< Tkey, TGeoMatrix * > transformation_matrices_aligned
Definition: PndLmdDim.h:664
bool operator<(const Tkey &comp) const
Definition: PndLmdDim.h:151
double end_seg_upstream
Definition: PndLmdDim.h:500
void Propagate(TVector3 &pos, TVector3 &momdir)
Definition: PndLmdDim.h:125
double plane_half_tilt_phi
Definition: PndLmdDim.h:307
double cvd_disc_rad
Definition: PndLmdDim.h:322
double phi_bend
Definition: PndLmdDim.h:504
double maps_passive_right
Definition: PndLmdDim.h:421
double maps_die_height
Definition: PndLmdDim.h:436
double * plane_pos_z
Definition: PndLmdDim.h:285
double cvd_tilt_theta
Definition: PndLmdDim.h:354
STL namespace.
double pipe_thickness
Definition: PndLmdDim.h:491
map< Tkey, TGeoMatrix * > transformation_matrices
Definition: PndLmdDim.h:663
double rad_entrance
Definition: PndLmdDim.h:485
double length_pipe
Definition: PndLmdDim.h:495
double sensor_offset_x
Definition: PndLmdDim.h:460
bool operator==(const Tkey &comp) const
Definition: PndLmdDim.h:180
double sensor_tilt_phi
Definition: PndLmdDim.h:466
double plane_half_tilt_psi
Definition: PndLmdDim.h:311
double side_offset_x
Definition: PndLmdDim.h:363
double outer_rad
Definition: PndLmdDim.h:318
double box_size_y
Definition: PndLmdDim.h:477
unsigned int i
Definition: P4_F32vec4.h:33
double side_tilt_psi
Definition: PndLmdDim.h:373
int Get_sensor_id(int ihalf, int iplane, int imodule, int iside, int idie, int isensor)
Definition: PndLmdDim.h:538
Tkey(string key)
Definition: PndLmdDim.h:205
double maps_passive_bottom
Definition: PndLmdDim.h:419
unsigned int sensIDoffset
Definition: PndLmdDim.h:269
double pol_side_dist_min
Definition: PndLmdDim.h:334
double length_transision
Definition: PndLmdDim.h:493
double sensor_offset_y
Definition: PndLmdDim.h:462
double pos_plane_0
Definition: PndLmdDim.h:497
double maps_width
Definition: PndLmdDim.h:427
double side_offset_y
Definition: PndLmdDim.h:365
signed char die
Definition: PndLmdDim.h:149
double cvd_disc_thick_half
Definition: PndLmdDim.h:324
Alignment Parameter Class for LMD.
double maps_active_offset_x
Definition: PndLmdDim.h:430
void Get_sensor_by_id(const int sensor_id, int &ihalf, int &iplane, int &imodule, int &iside, int &idie, int &isensor)
Definition: PndLmdDim.h:547
double kapton_disc_thick_half
Definition: PndLmdDim.h:360
double rot_y
Definition: PndLmdDim.h:515
double rad_pipe
Definition: PndLmdDim.h:489
double end_seg_bend
Definition: PndLmdDim.h:509
double box_thickness
Definition: PndLmdDim.h:481
double pol_side_lg_half
Definition: PndLmdDim.h:332
double pos_rib
Definition: PndLmdDim.h:483
double half_tilt_theta
Definition: PndLmdDim.h:296
double cvd_offset_x
Definition: PndLmdDim.h:346
double maps_passive_left
Definition: PndLmdDim.h:420
double maps_active_offset_y
Definition: PndLmdDim.h:431
double cvd_tilt_phi
Definition: PndLmdDim.h:352
Tkey()
Definition: PndLmdDim.h:235
double box_size_z
Definition: PndLmdDim.h:479
double die_offset_x
Definition: PndLmdDim.h:447
signed char sensor
Definition: PndLmdDim.h:150
double cvd_disc_even_odd_offset
Definition: PndLmdDim.h:327
int maps_n_col
Definition: PndLmdDim.h:409
double pos_rot_z
Definition: PndLmdDim.h:506
double side_tilt_theta
Definition: PndLmdDim.h:371
double cvd_disc_dist
Definition: PndLmdDim.h:338
double die_gap
Definition: PndLmdDim.h:433
map< Tkey, vector< double > > offsets
Definition: PndLmdDim.h:580
void Set(const double *matrix)
Definition: PndLmdDim.h:88
int maps_n_row
Definition: PndLmdDim.h:410
bool ** enabled
Definition: PndLmdDim.h:412
double cvd_offset_z
Definition: PndLmdDim.h:350
double inner_rad
Definition: PndLmdDim.h:316
Tkey(const Tkey &copy)
Definition: PndLmdDim.h:185
double maps_die_width
Definition: PndLmdDim.h:435
double plane_half_offset_x
Definition: PndLmdDim.h:301
PndLmdDimPropMat(double mom_init)
Definition: PndLmdDim.h:97
double maps_active_height
Definition: PndLmdDim.h:423
int n_cvd_discs
Definition: PndLmdDim.h:320
signed char module
Definition: PndLmdDim.h:147
double pi
Definition: PndLmdDim.h:279
int sign(T val)
Definition: PndCADef.h:61
unsigned int nmodules
Definition: PndLmdDim.h:283
double half_offset_x
Definition: PndLmdDim.h:288
double plane_half_offset_z
Definition: PndLmdDim.h:305
double half_tilt_psi
Definition: PndLmdDim.h:298
double half_offset_y
Definition: PndLmdDim.h:290
double rot_x
Definition: PndLmdDim.h:514
double maps_height
Definition: PndLmdDim.h:428
int n_sensors
Definition: PndLmdDim.h:414
signed char plane
Definition: PndLmdDim.h:146
double cvd_tilt_psi
Definition: PndLmdDim.h:356
double side_tilt_phi
Definition: PndLmdDim.h:369
double sensor_offset_z
Definition: PndLmdDim.h:464
double sensor_tilt_theta
Definition: PndLmdDim.h:468
double rad_exit
Definition: PndLmdDim.h:487
double maps_thickness
Definition: PndLmdDim.h:417
map< double, PndLmdDimPropMat > propagation_matrices
Definition: PndLmdDim.h:775
signed char side
Definition: PndLmdDim.h:148
int Generate_keynumber(unsigned int ihalf=0, unsigned int iplane=0, unsigned int imodule=0, unsigned int iside=0, unsigned int idie=0, unsigned int isensor=0)
Definition: PndLmdDim.h:642
double maps_active_width
Definition: PndLmdDim.h:422
unsigned int n_planes
Definition: PndLmdDim.h:281
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:64
double pos_x
Definition: PndLmdDim.h:511
double box_size_x
Definition: PndLmdDim.h:475
double maps_passive_top
Definition: PndLmdDim.h:418
signed char half
Definition: PndLmdDim.h:145
double side_offset_z
Definition: PndLmdDim.h:367
double die_offset_z
Definition: PndLmdDim.h:451