FairRoot/PandaRoot
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
PndEmcStructure Class Reference

geometry helper class More...

#include <PndEmcStructure.h>

Inheritance diagram for PndEmcStructure:

Public Types

typedef std::map< Int_t, Float_t > mapper
 

Public Member Functions

virtual ~PndEmcStructure ()
 
const mapperGetEmcX () const
 
const mapperGetEmcY () const
 
const mapperGetEmcZ () const
 
const PndEmcTciXtalMapGetTciXtalMap () const
 
PndEmcTwoCoordIndexlocateIndex (double theta, double phi) const
 
void Print (string, Int_t option=1) const
 

Static Public Member Functions

static PndEmcStructureInstance ()
 
static PndEmcStructureInstance (TGeoManager *)
 

Protected Member Functions

 PndEmcStructure (TGeoManager *)
 

Private Member Functions

bool crystal_name_analysis (TString, int &module, int &copy, int &row, int &crystal)
 

Private Attributes

mapper emcX
 
mapper emcY
 
mapper emcZ
 
PndEmcTciXtalMap fTciXtalMap
 

Static Private Attributes

static PndEmcStructure_instance = 0
 

Detailed Description

geometry helper class

Definition at line 25 of file PndEmcStructure.h.

Member Typedef Documentation

typedef std::map<Int_t, Float_t> PndEmcStructure::mapper

Definition at line 29 of file PndEmcStructure.h.

Constructor & Destructor Documentation

PndEmcStructure::~PndEmcStructure ( )
virtual

Definition at line 52 of file PndEmcStructure.cxx.

52 {}
PndEmcStructure::PndEmcStructure ( TGeoManager *  geoMan)
protected

Definition at line 84 of file PndEmcStructure.cxx.

References CAMath::ATan2(), crystal_name_analysis(), Double_t, dz, emcX, emcY, emcZ, fTciXtalMap, PndEmcMapper::GetTCI(), h1, h2, PndEmcMapper::Instance(), next, pos, rot, row, CAMath::Sqrt(), trans, and TString.

84  : emcX(), emcY(), emcZ(), fTciXtalMap()
85 {
86 
87  // Instantiate mapper to convert from detId to TCI
88  // PndEmcMapper should be initiated before
90 
91  int module,copy,row,crystal;
92  TString node_path;
93  //const char* crystal_name; //[R.K. 01/2017] unused variable?
94  const TGeoMatrix *crystal_matrix;
95  const double *trans;
96  TGeoRotation geoRot;
97  TRotation rot;
98  //Double_t phi, theta, psi; //[R.K. 01/2017] unused variable?
99  int detId;
100 
101  TGeoIterator next(geoMan->GetTopVolume());
102  TGeoNode *node;
103 
104  while ((node=next())) {
105  next.GetPath(node_path);
106 
107  bool isEmcModule=crystal_name_analysis(node_path,module,copy,row,crystal);
108  if (!isEmcModule) continue;
109 
110 
111  detId = module*100000000 + row*1000000 + copy*10000 + crystal;
112 
113 
114  PndEmcTwoCoordIndex *tci=fEmcMap->GetTCI(detId);
115 
116  if (tci==0)
117  {
118  cout<<"Not found tci for index = "<<detId<<" in PndEmcStructure"<<endl;
119  cout<<node_path<<endl;
120  abort();
121  }
122 
123  // Obtaine TGeoMatrix matrix for the current node, later it is factorised to translation and rotation
124  crystal_matrix=next.GetCurrentMatrix();
125 
126  trans=crystal_matrix->GetTranslation();
127  TVector3 pos(trans[0],trans[1],trans[2]);
128  emcX[detId] = pos.X();
129  emcY[detId] = pos.Y();
130  emcZ[detId] = pos.Z();
131 
132  geoRot.SetMatrix(crystal_matrix->GetRotationMatrix());
133 
134  PndEmcXtal *xtal;
135 
136  // Check the geometry type and convert to TGeoTrap.
137 
138  TString shapeType = node->GetVolume()->GetShape()->ClassName();
139  if (shapeType == "TGeoTrap") {
140  TGeoTrap *trap = dynamic_cast<TGeoTrap *>(
141  node->GetVolume()->GetShape());
142 
143  xtal = new PndEmcXtal(tci, *trap, pos, geoRot);
144  } else if (shapeType == "TGeoArb8") {
145  TGeoArb8 *arb8 = dynamic_cast<TGeoArb8 *>(
146  node->GetVolume()->GetShape());
147 
148  // Approximate the shape with a TGeoTrap.
149 
150  Double_t *verts = arb8->GetVertices();
151 
152  Double_t dz = arb8->GetDz();
153  Double_t tx = (verts[4*2+0] + verts[5*2+0] + verts[6*2+0] + verts[7*2+0] - verts[0*2+0] - verts[1*2+0] - verts[2*2+0] - verts[3*2+0]) / (2 * 4 * dz);
154  Double_t ty = (verts[4*2+1] + verts[5*2+1] + verts[6*2+1] + verts[7*2+1] - verts[0*2+1] - verts[1*2+1] - verts[2*2+1] - verts[3*2+1]) / (2 * 4 * dz);
155 
156  Double_t thetac = TMath::ATan(TMath::Sqrt(tx*tx + ty*ty)) * TMath::RadToDeg();
157  Double_t phic = TMath::ATan2(ty, tx) * TMath::RadToDeg();
158 
159  Double_t h1 = (verts[1*2+1] + verts[2*2+1] - verts[0*2+1] - verts[3*2+1])/(2*2);
160  Double_t bl1 = (verts[3*2+0] - verts[0*2+0]) / 2;
161  Double_t tl1 = (verts[2*2+0] - verts[1*2+0]) / 2;
162  Double_t alpha1 = TMath::ATan((verts[1*2+0]+verts[2*2+0]-verts[0*2+0]-verts[3*2+0])/(2*2*h1)) * TMath::RadToDeg();
163 
164  Double_t h2 = (verts[5*2+1] + verts[6*2+1] - verts[4*2+1] - verts[7*2+1])/(2*2);
165  Double_t bl2 = (verts[7*2+0] - verts[4*2+0]) / 2;
166  Double_t tl2 = (verts[6*2+0] - verts[5*2+0]) / 2;
167  Double_t alpha2 = TMath::ATan((verts[5*2+0]+verts[6*2+0]-verts[4*2+0]-verts[7*2+0])/(2*2*h2)) * TMath::RadToDeg();
168 
169  TGeoTrap *crystal_shape = new TGeoTrap(dz, thetac, phic,
170  h1, bl1, tl1, alpha1,
171  h2, bl2, tl2, alpha2);
172 
173  // Check the conversion.
174  //Double_t *verts2 = crystal_shape.GetVertices();
175  //for (size_t i = 0; i < 8; ++i) {
176  // if (TMath::Abs(verts[i*2+0] - verts2[i*2+0]) > 1e-8
177  // || TMath::Abs(verts[i*2+1] - verts2[i*2+1]) > 1e-8) {
178  // cout << "TGeoArb8 conversion bad in module " << module << " " << " x " << verts[i*2+0] << "->" << verts2[i*2+0] << " y " << verts[i*2+1] << "->" << verts2[i*2+1] << endl;
179  // }
180  //}
181 
182  xtal = new PndEmcXtal(tci, *crystal_shape, pos, geoRot);
183  } else if (shapeType == "TGeoBBox" || shapeType == "TGeoScaledShape") {
184  TGeoBBox const *box = dynamic_cast<TGeoBBox const*>(
185  node->GetVolume()->GetShape());
186 
187  // Convert to TGeoTrap.
188  TGeoTrap *crystal_shape = new TGeoTrap(box->GetDZ(), 0, 0,
189  box->GetDY(), box->GetDX(), box->GetDX(), 0,
190  box->GetDY(), box->GetDX(), box->GetDX(), 0);
191 
192  xtal = new PndEmcXtal(tci,*crystal_shape,pos,geoRot);
193 
194  } else {
195  cout << "Unknown geometry type " << shapeType << " in module " << module << endl;
196  abort();
197  }
198 
199  fTciXtalMap[tci]=xtal;
200  }
201 }
int row
Definition: anaLmdDigi.C:67
TVector3 pos
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TGeoManager * geoMan
PndEmcTciXtalMap fTciXtalMap
represents coordinates of one crystal
Definition: PndEmcXtal.h:36
stores crystal index coordinates (x,y) or (theta,phi)
PndEmcTwoCoordIndex * GetTCI(Int_t DetectorId)
Emc geometry mapper.
Definition: PndEmcMapper.h:22
TGeoTranslation * trans
bool crystal_name_analysis(TString, int &module, int &copy, int &row, int &crystal)
Double_t
static T ATan2(const T &y, const T &x)
TGeoRotation rot
static int next[96]
Definition: ranlxd.cxx:374
static PndEmcMapper * Instance()

Member Function Documentation

bool PndEmcStructure::crystal_name_analysis ( TString  node_path,
int &  module,
int &  copy,
int &  row,
int &  crystal 
)
private

Definition at line 203 of file PndEmcStructure.cxx.

References col, fabs(), getsign(), next, and row.

Referenced by PndEmcStructure().

204 {
205  // The following code extracts information about module, copy, row and crystal number from the path name of the TGeoNode which corresponds to emc crystal
206  // Name convention is according to PndEmc.cxx
207 
208 
209  if (!(node_path.Contains("emc") || node_path.Contains("Crystal") || node_path.Contains("FscModuleVolume"))) return false;
210  //cout << "node_path= " << node_path << endl;
212  // Case of old Fsc
214  if (node_path.Contains("Fsc_")){
215  //at the moment all the layers of module in Fsc are not taken into account, only the whole block is selected
216  if (node_path.Contains("FscLayer"))
217  return false;
218  sscanf(node_path.Data(),"cave/Fsc_%d/emc%dr%dc%d_0",&copy,&module,&row,&crystal);
219  return true;
220  }
221 
223  // Case of new Fsc
225  if (node_path.Contains("FscModuleVolume")){
226  //at the moment all the layers of module in Fsc are not taken into account, only the whole block is selected
227 // if (node_path.Contains("FscLayer") || node_path.Contains("FscTyvek") || node_path.Contains("FscFibHole"))
228  if (node_path.Contains("FscLayer") || node_path.Contains("FscFibHole"))
229  return false;
230 // cout << "node_path= " << node_path << endl;
231  int SupModCopy=0;
232  int LocCopy=0;
233  int dummyTyv = 0;
234  Int_t nSupCol=0;
235  Int_t nSupRow = 0;
236  Int_t nModCol = 0;
237  Int_t nModRow = 0;
238 
239 // sscanf(node_path.Data(),"cave/Emc%d_%d/FscModuleVolume_%d",&module,&copy,&ModCopy);
240  sscanf(node_path.Data(),"cave/Emc%d_%d/FscSuperModuleVolume_%d/FscTyvekVolume_%d/FscModuleVolume_%d",
241  &module,&copy,&SupModCopy,&dummyTyv,&LocCopy);
242  copy+=1;
243 
244  nSupCol = SupModCopy/100;
245  nSupRow = SupModCopy%100;
246 
247 
248  nModCol = LocCopy%2;
249  nModRow = LocCopy/2;
250 
251  crystal = (nSupCol - 1)*2 + nModCol + 1;
252  row = (nSupRow - 1)*2 + nModRow + 1;
253 
254 
255 
256  return true;
257  }
258 
260  // Case of barrel section with hole
262  if (node_path.Contains("Emc12Hole")){
263  //int tmp1,tmp2; //[R.K. 01/2017] unused variable?
264  if (node_path.Contains("EmcLayer2Hole"))
265  {
266  sscanf(node_path.Data(),"cave/Emc12Hole_%d/EmcLayer2Hole_0/emc%dr%dc%d_0$",&copy,&module,&row,&crystal);
267  }
268  else
269  {
270  sscanf(node_path.Data(),"cave/Emc12Hole_%d/EmcLayer1_0/emc%dr%dc%d_0$",&copy,&module,&row,&crystal);
271  }
272  return true;
273  }
274 
276  // Case of barrel
278  if (node_path.Contains("EmcLayer")){
279  int tmp1,tmp2;
280  sscanf(node_path.Data(),"cave/Emc%d_%d/EmcLayer%d_0/emc%dr%dc%d_0$",&tmp1,&copy,&tmp2,&module,&row,&crystal);
281  }
282 
284  // Case of new version of barrel section with hole
286  if (node_path.Contains("Slice_target")){
287  int iSM, copySM, iMod, copyMod, iAlv, iCry, copyCry;
288  char sgnMod, sgnCry, typeCry;
289  sscanf(node_path.Data(),
290  "cave/BarrelEMC_0/Slice_target_%d/SuperModule%d_Target_%d/Module%d%c_%d/Crystal-%d%c-%c%d_%d",
291  &copy, &iSM, &copySM, &iMod, &sgnMod, &copyMod, &iAlv, &sgnCry, &typeCry, &iCry, &copyCry);
292 
293  // nMod
294  if (sgnCry == 'p') module = 1;
295  else if (sgnCry == 'm') module = 2;
296  else return false;
297 
298  // nRow
299  row = (iAlv - 1) * 4 + iCry;
300 
301  // nCrys
302  if (module == 1) {
303  if (typeCry == 'L') crystal = (copyCry - 1) * 2 + 1;
304  if (typeCry == 'R') crystal = (copyCry - 1) * 2 + 2;
305  }
306  else if (module == 2) {
307  if (typeCry == 'R') crystal = (copyCry - 1) * 2 + 1;
308  if (typeCry == 'L') crystal = (copyCry - 1) * 2 + 2;
309  }
310  else {
311  std::cout << "Error!!!" << std::endl;
312  return false;
313  }
314  return true;
315  }
316 
318  // Case of new version of barrel
320  if (node_path.Contains("Slice")){
321  int iSM, copySM, iMod, copyMod, iAlv, iCry, copyCry;
322  char sgnMod, sgnCry, typeCry;
323  sscanf(node_path.Data(),
324  "cave/BarrelEMC_0/Slice_%d/SuperModule%d_%d/Module%d%c_%d/Crystal-%d%c-%c%d_%d",
325  &copy, &iSM, &copySM, &iMod, &sgnMod, &copyMod, &iAlv, &sgnCry, &typeCry, &iCry, &copyCry);
326 
327  // nMod
328  if (sgnCry == 'p') module = 1;
329  else if (sgnCry == 'm') module = 2;
330  else return false;
331 
332  // nRow
333  row = (iAlv - 1) * 4 + iCry;
334 
335  // nCrys
336  if (module == 1) {
337  if (typeCry == 'L') crystal = (copyCry - 1) * 2 + 1;
338  if (typeCry == 'R') crystal = (copyCry - 1) * 2 + 2;
339  }
340  else if (module == 2) {
341  if (typeCry == 'R') crystal = (copyCry - 1) * 2 + 1;
342  if (typeCry == 'L') crystal = (copyCry - 1) * 2 + 2;
343  }
344  else {
345  std::cout << "Error!!!" << std::endl;
346  return false;
347  }
348  return true;
349  }
350 
352  //case of test calorimeter
354  else if (node_path.Contains("EmcTest")) {
355  sscanf(node_path.Data(),"cave/EmcTest_%d/emc%dr%dc%d_0$",&copy,&module,&row,&crystal);
356  }
357 
359  //case new version of forward end-cap (from "emc_module3_2011_new.root" file)
361  else if (node_path.Contains("SubunitVolFwEndCap")) {
362  Int_t copyNoSub, copyNoBox, copyNoCrys;
363 
364  if (node_path.Contains("HalfSubunitVolFwEndCap")) //reading the HalfSubunits of FwEndCap
365  sscanf(node_path.Data(),"cave/Emc3_0/HalfSubunitVolFwEndCap_%d/BoxVol_%d/CrystalVol_%d", &copyNoSub, &copyNoBox, &copyNoCrys);
366  else //reading the Subunits of FwEndCap
367  sscanf(node_path.Data(),"cave/Emc3_0/SubunitVolFwEndCap_%d/BoxVol_%d/CrystalVol_%d", &copyNoSub, &copyNoBox, &copyNoCrys);
368 
369  //cout << "copyNoSub = "<< copyNoSub<<" ,copyNoBox = "<< copyNoBox<<" ,copyNoCrys = "<< copyNoCrys<<"\n";
370 
371  //some presetting of parameters:
372  Int_t SubunitCol=-100, SubunitRow=-100, CrystalCol=-100, CrystalRow=-100;
373  Int_t RestOfHalfSubunitRowNo[6] = {5, 6, 7, 8, 9, 9};//row number of 6 peripheral half Subunits
374  Int_t RestOfHalfSubunitColNo[6] = {8, 7, 6, 5, 3, 2};//column number of 6 peripheral half Subunits
375  //determination of SubunitRow and SubunitCol:
376  if (copyNoSub >= 1 && copyNoSub <= 10){ // the middle row of 10 HalfSubunits
377  SubunitRow = 0;
378  SubunitCol = pow(-1.,1+(copyNoSub-1)/5)*(4+(copyNoSub-1)%5);
379  }
380  else if (copyNoSub > 11 && copyNoSub < 24){ // the upper and lower pipe-region row of 10 HalfSubunits
381  SubunitRow = pow(-1.,(copyNoSub-11)/7)*2;
382  SubunitCol = ((copyNoSub-11)%7)-3;
383  }
384  else if ((copyNoSub >= 25 && copyNoSub <= 27) || copyNoSub==57 || copyNoSub==58){ // the outermost left column of 5 HalfSubunits
385  if (copyNoSub >= 25 && copyNoSub <= 27)
386  SubunitRow = copyNoSub - 25;
387  else
388  SubunitRow = copyNoSub - 59;
389  SubunitCol = -9;
390  }
391  else if (copyNoSub >= 40 && copyNoSub <= 44){ // the outermost right column of 5 HalfSubunits
392  SubunitRow = 42 - copyNoSub;
393  SubunitCol = 9;
394  }
395  else if (copyNoSub == 28 || copyNoSub ==39 || copyNoSub == 45 || copyNoSub == 56) { // 4 single peripheral half Subunits
396  if (copyNoSub <= 39)
397  SubunitRow = RestOfHalfSubunitRowNo[0];
398  else
399  SubunitRow = -RestOfHalfSubunitRowNo[0];
400  SubunitCol = pow(-1.,copyNoSub-1)*RestOfHalfSubunitColNo[0];
401  }
402  else if (copyNoSub == 29 || copyNoSub == 38 || copyNoSub == 46 || copyNoSub == 55) { // 4 single peripheral half Subunits
403  if (copyNoSub <= 38)
404  SubunitRow = RestOfHalfSubunitRowNo[1];
405  else
406  SubunitRow = -RestOfHalfSubunitRowNo[1];
407  SubunitCol = pow(-1.,copyNoSub)*RestOfHalfSubunitColNo[1];
408  }
409  else if (copyNoSub == 30 || copyNoSub == 37 || copyNoSub == 47 || copyNoSub == 54) { // 4 single peripheral half Subunits
410  if (copyNoSub <= 37)
411  SubunitRow = RestOfHalfSubunitRowNo[2];
412  else
413  SubunitRow = -RestOfHalfSubunitRowNo[2];
414  SubunitCol = pow(-1.,copyNoSub-1)*RestOfHalfSubunitColNo[2];
415  }
416  else if (copyNoSub == 31 || copyNoSub == 36 || copyNoSub == 48 || copyNoSub == 53) { // 4 single peripheral half Subunits
417  if (copyNoSub <= 36)
418  SubunitRow = RestOfHalfSubunitRowNo[3];
419  else
420  SubunitRow = -RestOfHalfSubunitRowNo[3];
421  SubunitCol = pow(-1.,copyNoSub)*RestOfHalfSubunitColNo[3];
422  }
423  else if (copyNoSub == 32 || copyNoSub == 35 || copyNoSub == 49 || copyNoSub == 52) { // 4 single peripheral half Subunits
424  if (copyNoSub <= 35)
425  SubunitRow = RestOfHalfSubunitRowNo[4];
426  else
427  SubunitRow = -RestOfHalfSubunitRowNo[4];
428  SubunitCol = pow(-1.,copyNoSub-1)*RestOfHalfSubunitColNo[4];
429  }
430  else if (copyNoSub == 33 || copyNoSub == 34 || copyNoSub == 50 || copyNoSub == 51) { // 4 single peripheral half Subunits
431  if (copyNoSub <= 34)
432  SubunitRow = RestOfHalfSubunitRowNo[5];
433  else
434  SubunitRow = -RestOfHalfSubunitRowNo[5];
435  SubunitCol = pow(-1.,copyNoSub)*RestOfHalfSubunitColNo[5];
436  }
437  else if (copyNoSub >= 61 && copyNoSub <= 74){ // the middle column of full Subunits
438  SubunitRow = pow(-1.,(copyNoSub-61)/7)*((copyNoSub - 61)%7 + 3);
439  SubunitCol = 0;
440  }
441  else if (copyNoSub%100 >= 1 && copyNoSub%100 <= 5){ // rest of the full Subunits
442  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
443  SubunitRow = 1;
444  SubunitCol = pow(-1.,1+copyNoSub/100)*(3 + copyNoSub%100);
445  }
446  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
447  SubunitRow = -1;
448  SubunitCol = pow(-1.,copyNoSub/100)*(3 + copyNoSub%100);
449  }
450  }
451  else if (copyNoSub%100 >= 6 && copyNoSub%100 <= 11){ // rest of the full Subunits
452  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
453  SubunitRow = 2;
454  SubunitCol = pow(-1.,1+copyNoSub/100)*(2 + copyNoSub%100 - 5);
455  }
456  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
457  SubunitRow = -2;
458  SubunitCol = pow(-1.,copyNoSub/100)*(2 + copyNoSub%100 - 5);
459  }
460  }
461  else if (copyNoSub%100 >= 12 && copyNoSub%100 <= 19){ // rest of the full Subunits
462  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
463  SubunitRow = 3;
464  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 11);
465  }
466  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
467  SubunitRow = -3;
468  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 11);
469  }
470  }
471  else if (copyNoSub%100 >= 20 && copyNoSub%100 <= 27){ // rest of the full Subunits
472  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
473  SubunitRow = 4;
474  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 19);
475  }
476  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
477  SubunitRow = -4;
478  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 19);
479  }
480  }
481  else if (copyNoSub%100 >= 28 && copyNoSub%100 <= 34){ // rest of the full Subunits
482  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
483  SubunitRow = 5;
484  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 27);
485  }
486  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
487  SubunitRow = -5;
488  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 27);
489  }
490  }
491  else if (copyNoSub%100 >= 35 && copyNoSub%100 <= 40){ // rest of the full Subunits
492  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
493  SubunitRow = 6;
494  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 34);
495  }
496  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
497  SubunitRow = -6;
498  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 34);
499  }
500  }
501  else if (copyNoSub%100 >= 41 && copyNoSub%100 <= 45){ // rest of the full Subunits
502  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
503  SubunitRow = 7;
504  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 40);
505  }
506  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
507  SubunitRow = -7;
508  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 40);
509  }
510  }
511  else if (copyNoSub%100 >= 46 && copyNoSub%100 <= 49){ // rest of the full Subunits
512  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
513  SubunitRow = 8;
514  SubunitCol = pow(-1.,1+copyNoSub/100)*(copyNoSub%100 - 45);
515  }
516  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
517  SubunitRow = -8;
518  SubunitCol = pow(-1.,copyNoSub/100)*(copyNoSub%100 - 45);
519  }
520  }
521  else if (copyNoSub%100 == 50 ){ // rest of the full Subunits
522  if (copyNoSub/100 == 1 || copyNoSub/100 == 4){
523  SubunitRow = 9;
524  SubunitCol = pow(-1.,1+copyNoSub/100);
525  }
526  else {//here necessarily: copyNoSub/100 == 2 || copyNoSub/100 == 3
527  SubunitRow = -9;
528  SubunitCol = pow(-1.,copyNoSub/100);
529  }
530  }
531 
532  Int_t flag=1; //this flag introducing is not necessary in this class but is also harmless
533  if (fabs(SubunitRow)<=1 && fabs(SubunitCol)<=3) flag=0; // empty copyNoSub in the beam-pipe area
534  else if (fabs(SubunitCol)==9 && fabs(SubunitRow)>=3) flag=0; // empty copyNoSub in the residual area
535  else if (fabs(SubunitCol)==8 && fabs(SubunitRow)>=6) flag=0; // -||-
536  else if (fabs(SubunitCol)==7 && fabs(SubunitRow)>=7) flag=0; // -||-
537  else if (fabs(SubunitCol)==6 && fabs(SubunitRow)>=8) flag=0; // -||-
538  else if (fabs(SubunitCol)>=4 && fabs(SubunitRow)==9) flag=0; // -||-
539 
540  //determination of CrystalRow and CrystalCol:
541  if (flag && copyNoSub >= 61){//determination of CrystalRow and CrystalCol for the 214 full Subunits
542  if(copyNoBox==1 || copyNoBox==4){
543  if(copyNoCrys==1 || copyNoCrys==3)
544  CrystalRow = SubunitRow*4 + (3-getsign(SubunitRow))/2;
545  else //here necessarily: copyNoCrys==2 || copyNoCrys==4
546  CrystalRow = SubunitRow*4 + (1-getsign(SubunitRow))/2;
547  }
548  else if(copyNoBox==2 || copyNoBox==3){
549  if(copyNoCrys==1 || copyNoCrys==3)
550  CrystalRow = SubunitRow*4 - (1+getsign(SubunitRow))/2;
551  else //here necessarily: copyNoCrys==2 || copyNoCrys==4
552  CrystalRow = SubunitRow*4 - (3+getsign(SubunitRow))/2;
553  }
554 
555  if (copyNoSub >= 61 && copyNoSub <= 74) {//determination of CrystalCol for the 14 middle column full Subunits
556  if(copyNoBox==1 || copyNoBox==3){
557  if(copyNoCrys==1 || copyNoCrys==4)
558  CrystalCol = 2;
559  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
560  CrystalCol = 1;
561  }
562  else if(copyNoBox==4 || copyNoBox==2){
563  if(copyNoCrys==1 || copyNoCrys==4)
564  CrystalCol = -1;
565  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
566  CrystalCol = -2;
567  }
568  }
569  else { //determination of CrystalCol for the other 4*50=200 full Subunits
570  if(copyNoBox==1 || copyNoBox==3){
571  if(copyNoCrys==1 || copyNoCrys==4)
572  CrystalCol = SubunitCol*4 + (3+getsign(SubunitCol))/2;
573  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
574  CrystalCol = SubunitCol*4 + (1+getsign(SubunitCol))/2;
575  }
576  else if(copyNoBox==4 || copyNoBox==2){
577  if(copyNoCrys==1 || copyNoCrys==4)
578  CrystalCol = SubunitCol*4 - (1-getsign(SubunitCol))/2;
579  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
580  CrystalCol = SubunitCol*4 - (3-getsign(SubunitCol))/2;
581  }
582  }
583  }
584 
585  else if (flag && copyNoSub <= 10) {//determination of CrystalRow and CrystalCol for the 10 middle-row half Subunits
586  if(copyNoCrys==1 || copyNoCrys==3)
587  CrystalRow = 1;
588  else //here necessarily: copyNoCrys==2 || copyNoCrys==4
589  CrystalRow = -1;
590 
591  if(copyNoBox==1){
592  if(copyNoCrys==1 || copyNoCrys==4)
593  CrystalCol = SubunitCol*4 + (3+getsign(SubunitCol))/2;
594  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
595  CrystalCol = SubunitCol*4 + (1+getsign(SubunitCol))/2;
596  }
597  else { //here necessarily: copyNoBox==2
598  if(copyNoCrys==1 || copyNoCrys==4)
599  CrystalCol = SubunitCol*4 - (1-getsign(SubunitCol))/2;
600  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
601  CrystalCol = SubunitCol*4 - (3-getsign(SubunitCol))/2;
602  }
603  }
604  else if (flag && copyNoSub > 11 && copyNoSub < 24) {//determination of CrystalRow and CrystalCol for the 10 upper and lower pipe-region row of half subunits
605  if(copyNoCrys==1 || copyNoCrys==3)
606  CrystalRow = (1+17*getsign(SubunitRow))/2;
607  else //here necessarily: copyNoCrys==2 || copyNoCrys==4
608  CrystalRow = (17*getsign(SubunitRow)-1)/2;
609 
610  if (SubunitCol == 0 && copyNoBox==1) {
611  if(copyNoCrys==1 || copyNoCrys==4)
612  CrystalCol = 2;
613  else
614  CrystalCol = 1;
615  }
616  else if (SubunitCol == 0) {//here necessarily: copyNoBox==2
617  if (copyNoCrys==1 || copyNoCrys==4)
618  CrystalCol = -1;
619  else
620  CrystalCol = -2;
621  }
622  else if(copyNoBox==1){
623  if(copyNoCrys==1 || copyNoCrys==4)
624  CrystalCol = SubunitCol*4 + (3+getsign(SubunitCol))/2;
625  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
626  CrystalCol = SubunitCol*4 + (1+getsign(SubunitCol))/2;
627  }
628  else { //here necessarily: copyNoBox==2
629  if(copyNoCrys==1 || copyNoCrys==4)
630  CrystalCol = SubunitCol*4 - (1-getsign(SubunitCol))/2;
631  else //here necessarily: copyNoCrys==3 || copyNoCrys==2
632  CrystalCol = SubunitCol*4 - (3-getsign(SubunitCol))/2;
633  }
634  }
635  else if (flag && ((copyNoSub >= 25 && copyNoSub <= 27) || copyNoSub == 57 || copyNoSub == 58 || (copyNoSub >= 40 && copyNoSub <= 44))) {//determination of CrystalRow and CrystalCol for the 2 outermost left and right columns of half Subunits
636  if(copyNoBox==1 && (copyNoCrys==1 || copyNoCrys==4))
637  if (copyNoSub == 25 || copyNoSub == 42)
638  CrystalRow = 2;
639  else
640  CrystalRow = SubunitRow*4 + (3+getsign(SubunitRow))/2;
641  else if (copyNoBox==1)//here necessarily: copyNoCrys==3 || copyNoCrys==2
642  if (copyNoSub == 25 || copyNoSub == 42)
643  CrystalRow = 1;
644  else
645  CrystalRow = SubunitRow*4 + (1+getsign(SubunitRow))/2;
646 
647  else if (copyNoBox==2 && (copyNoCrys==1 || copyNoCrys==4))
648  if (copyNoSub == 25 || copyNoSub == 42)
649  CrystalRow = -1;
650  else
651  CrystalRow = SubunitRow*4 - (1-getsign(SubunitRow))/2;
652  else if (copyNoBox==2)//here necessarily: copyNoCrys==3 || copyNoCrys==2 // FIXME Use braces for explicit logic, possible ambiguity
653  if (copyNoSub == 25 || copyNoSub == 42)
654  CrystalRow = -2;
655  else
656  CrystalRow = SubunitRow*4 - (3-getsign(SubunitRow))/2;
657 
658  if (copyNoCrys==1 || copyNoCrys==3)
659  if (!(copyNoSub >= 40 && copyNoSub <= 44))
660  CrystalCol = -36;
661  else
662  CrystalCol = 35;
663  else //here necessarily: copyNoCrys==4 || copyNoCrys==2
664  if (!(copyNoSub >= 40 && copyNoSub <= 44))
665  CrystalCol = -35;
666  else
667  CrystalCol = 36;
668  }
669  else if (flag && (copyNoSub == 28 || copyNoSub ==29 || copyNoSub == 38 || copyNoSub == 39 || copyNoSub == 45 || copyNoSub == 46 || copyNoSub == 55 || copyNoSub == 56)) {//determination of CrystalRow and CrystalCol for the 8 single peripheral half Subunits
670  if (copyNoBox==1 && (copyNoCrys==1 || copyNoCrys==4))
671  CrystalRow = SubunitRow*4 + (3-getsign(SubunitRow))/2;
672  else if (copyNoBox==1) //here necessarily: copyNoCrys==3 || copyNoCrys==2
673  CrystalRow = SubunitRow*4 + (1-getsign(SubunitRow))/2;
674  else if(copyNoCrys==1 || copyNoCrys==4) //here necessarily: copyNoBox==2
675  CrystalRow = SubunitRow*4 - (1+getsign(SubunitRow))/2;
676  else //here necessarily: copyNoBox==2 && (copyNoCrys==3 || copyNoCrys==2)
677  CrystalRow = SubunitRow*4 - (3+getsign(SubunitRow))/2;
678 
679  if (copyNoCrys==4 || copyNoCrys==2)
680  CrystalCol = SubunitCol*4 + (1-getsign(SubunitCol))/2;
681  else //here necessarily: copyNoCrys==3 || copyNoCrys==1
682  CrystalCol = SubunitCol*4 - (1+getsign(SubunitCol))/2;
683  }
684  else if (flag && ((copyNoSub >= 30 && copyNoSub <= 37) || (copyNoSub >=47 && copyNoSub <= 54))) {//determination of CrystalRow and CrystalCol for the 16 single peripheral half Subunits
685  if (copyNoCrys==1 || copyNoCrys==3)
686  CrystalRow = SubunitRow*4 + (1-3*getsign(SubunitRow))/2;
687  else //here necessarily: copyNoCrys==4 || copyNoCrys==2
688  CrystalRow = SubunitRow*4 - (1+3*getsign(SubunitRow))/2;
689 
690  if (copyNoBox==1 && (copyNoCrys==1 || copyNoCrys==4))
691  CrystalCol = SubunitCol*4 + (3+getsign(SubunitCol))/2;
692  else if (copyNoBox==1) //here necessarily: copyNoCrys==3 || copyNoCrys==2
693  CrystalCol = SubunitCol*4 + (1+getsign(SubunitCol))/2;
694  else if (copyNoCrys==1 || copyNoCrys==4) //here necessarily: copyNoBox==2
695  CrystalCol = SubunitCol*4 - (1-getsign(SubunitCol))/2;
696  else //here necessarily: copyNoBox==2 && (copyNoCrys==3 || copyNoCrys==2)
697  CrystalCol = SubunitCol*4 - (3-getsign(SubunitCol))/2;
698  }
699 
700 
701  module = 3;
702  crystal = CrystalCol + 36;// to avoid negative values for "crystal", since it's gonna be defatorized from detID in an easy way (see PndEmcPoint.cxx)
703  row = CrystalRow + 37;// to avoid negative values for "row"
704  copy = 0;
705 
706  }
707 
708  /*else if (node_path.Contains("CrystalVol_block")) {
709  module = 3;
710  copy = 999;
711  crystal = 999;
712  row = 999;
713  }*/
714 
716  //case new version of backward end-cap (from "emc_module4new.root" file) - 9.10.2008
718  else if (node_path.Contains("Quarter4Vol")) {
719 
720  int copyNoSub, copyNoBox, copyNoCrys;
721  int tmp1, tmp2;
722 
723  if (node_path.Contains("SubunitVol_")){
724  sscanf(node_path.Data(),"cave/Emc4_0/Quarter4Vol_%d/SubunitVol_%d/BoxVol_%d/CrystalVol_%d",&copy, &copyNoSub, &copyNoBox, &copyNoCrys);
725  copyNoSub-=1;
726  }else{
727  if (node_path.Contains("BoxVol_")){
728  sscanf(node_path.Data(),"cave/Emc4_0/Quarter4Vol_%d/SubunitVol%d_%d/BoxVol_%d/CrystalVol_%d",&copy,&tmp1, &copyNoSub, &copyNoBox, &copyNoCrys);
729  copyNoSub-=1;
730  }else{
731  sscanf(node_path.Data(),"cave/Emc4_0/Quarter4Vol_%d/SubunitVol%d_%d/BoxVol%d_%d/CrystalVol_%d",&copy, &tmp1, &copyNoSub, &tmp2, &copyNoBox, &copyNoCrys);
732  copyNoSub-=1;
733  }
734  }
735 
736  Int_t col=0, nRow=-1, nCrys=-1; //k1=0, //[R.K. 01/2017] unused variable?
737  Int_t subrow=4; // 4 crystals in each subvolume
738  Int_t next=0; // starts (from the middle) next column
739 
740  // Below: for BwEndCap we have 13 subunits
741  // Now, 26.02.2009, 3 crystals in the middle's subunit are added =>
742  // => number of Subunits for BwEncCap & straight geometry is the same
743  if((copyNoSub >= 0) && (copyNoSub <= 3)){
744  next = copyNoSub;
745  col = 0;
746  }else if((copyNoSub >= 4) && (copyNoSub <= 7)){
747  next = (copyNoSub-4);
748  col = 1;
749  }else if((copyNoSub >= 8) && (copyNoSub <= 10)){
750  next = (copyNoSub-8);
751  col = 2;
752  }else if((copyNoSub >= 11) && (copyNoSub <= 12)){
753  next = (copyNoSub-11);
754  col = 3;
755  }
756 
757  Int_t flag4=1;
758  // 26.02.2009
759  // "next" means "row", "col" means "col"
760 
761  if (next>1 && col>2) flag4=0; // corner's subunit + one below
762  if (next>2 && col>1) flag4=0; // + one from the corner to the left
763 
764  if (flag4!=0){
765  if ( (copyNoBox == 0) || (copyNoBox == 3) ){
766  if(copyNoCrys == 1 || copyNoCrys == 3){
767  nCrys = next*4 + 3;
768  }else if (copyNoCrys == 0 || copyNoCrys == 2){
769  nCrys = next*4 + 4;
770  }
771  }else if ( (copyNoBox == 1) || (copyNoBox == 2) ){
772  if(copyNoCrys == 0 || copyNoCrys == 2){
773  nCrys = next*4 + 2;
774  }else if (copyNoCrys == 1 || copyNoCrys == 3){
775  nCrys = next*4 + 1;
776  }
777  }
778  if ( (copyNoBox == 0) || (copyNoBox == 2) ){
779  if(copyNoCrys == 0 || copyNoCrys == 3){
780  nRow = subrow*col + 4;
781  }else if (copyNoCrys == 1 || copyNoCrys == 2){
782  nRow = subrow*col + 3;
783  }
784  }else if ( (copyNoBox == 1) || (copyNoBox == 3) ){
785  if(copyNoCrys == 0 || copyNoCrys == 3){
786  nRow = subrow*col + 2;
787  }else if (copyNoCrys == 1 || copyNoCrys == 2){
788  nRow = subrow*col + 1;
789  }
790  }
791  }
792  module = 4;
793  row = nRow;
794  crystal = nCrys;
795 
796  }
797 
799  // backward endcap - emc_module4_StraightGeo24.4.root
801  else if (node_path.Contains("QuarterNewVol")){
802  module = 4;
803  int copyNoSub, copyNoBox, copyNoCrys;
804  int tmp1, tmp2;
805 
806  if (node_path.Contains("SubunitVol_")){
807  sscanf(node_path.Data(),"cave/Emc4_0/QuarterNewVol_%d/SubunitVol_%d/BoxVol_%d/CrystalVol_%d",&copy, &copyNoSub, &copyNoBox, &copyNoCrys);
808  copyNoSub-=1;
809  }else{
810  if (node_path.Contains("BoxVol_")){
811  sscanf(node_path.Data(),"cave/Emc4_0/QuarterNewVol_%d/SubunitVol%d_%d/BoxVol_%d/CrystalVol_%d",&copy, &tmp1, &copyNoSub, &copyNoBox, &copyNoCrys);
812  copyNoSub-=1;
813  }else{
814  sscanf(node_path.Data(),"cave/Emc4_0/QuarterNewVol_%d/SubunitVol%d_%d/BoxVol%d_%d/CrystalVol_%d",&copy, &tmp1, &copyNoSub, &tmp2, &copyNoBox, &copyNoCrys);
815  copyNoSub-=1;
816  }
817  }
818 
819  Int_t col=0, nRow=-1, nCrys=-1;
820  Int_t next=0; // starts (from the middle) next column
821 
822  // Below: for BwEndCap we have 13 subunits
823  // Now, 26.02.2009, 3 crystals in the middle's subunit are added =>
824  // => number of Subunits for BwEncCap & straight geometry is the same
825  if((copyNoSub >= 0) && (copyNoSub <= 2)){
826  next = copyNoSub+1;
827  col = 0;
828  }else if((copyNoSub >= 3) && (copyNoSub <= 6)){
829  next = (copyNoSub-3);
830  col = 1;
831  }else if((copyNoSub >= 7) && (copyNoSub <= 10)){
832  next = (copyNoSub-7);
833  col = 2;
834  }else if((copyNoSub >= 11) && (copyNoSub <= 13)){
835  next = (copyNoSub-11);
836  col = 3;
837  }
838 // cout << "copyNoSub= " << copyNoSub << " copyNoBoxs= " << copyNoBox << " copyNoCrys= " << copyNoCrys << endl;
839 
840  Int_t flag4=1;
841  // 26.02.2009
842  // "next" means "row", "col" means "col"
843 
844  if (next==3 && col==3) flag4=0; // corner's subunit + one below
845  if (next==0 && col==0) flag4=0; // + one from the corner to the left
846 // cout << "next= " << next << " col= " << col << endl;
847  if (flag4!=0){
848  if ( (copyNoBox == 0) || (copyNoBox == 3) ){
849  if(copyNoCrys == 1 || copyNoCrys == 3){
850  nCrys = next*4 + 3;
851  }else if (copyNoCrys == 0 || copyNoCrys == 2){
852  nCrys = next*4 + 4;
853  }
854  }else if ( (copyNoBox == 1) || (copyNoBox == 2) ){
855  if(copyNoCrys == 0 || copyNoCrys == 2){
856  nCrys = next*4 + 2;
857  }else if (copyNoCrys == 1 || copyNoCrys == 3){
858  nCrys = next*4 + 1;
859  }
860  }
861  if ( (copyNoBox == 0) || (copyNoBox == 2) ){
862  if(copyNoCrys == 0 || copyNoCrys == 3){
863  nRow = col*4 + 4;
864  }else if (copyNoCrys == 1 || copyNoCrys == 2){
865  nRow = col*4 + 3;
866  }
867  }else if ( (copyNoBox == 1) || (copyNoBox == 3) ){
868  if(copyNoCrys == 0 || copyNoCrys == 3){
869  nRow = col*4 + 2;
870  }else if (copyNoCrys == 1 || copyNoCrys == 2){
871  nRow = col*4 + 1;
872  }
873  }
874  }
875  module = 4;
876  row = nRow;
877  crystal = nCrys;
878  }
880  // Bwd Endcap geometry - 2017 version
881  // copy: quarter no. (0-4), row: submodule no. (0-9), crystal: cryst. no. (0-15)
883  else if(node_path.Contains("BWECinnerVol")){
884  int submodType = -1;
885  const char pathform[]="cave/Emc4_0/BWECouterVol_0/BWECinnerVol_0/BWECquarter_%d/BWECsubmodule%d_%d/PWOCrystal_%d";
886  sscanf(node_path.Data(),pathform, &copy, &submodType, &row, &crystal);
887  module = 4;
888  }
889 
892  else if(node_path.Contains("EmcProto")){
893  module=7;
894  copy=1;
895  sscanf(node_path.Data(),"cave/EmcProto_0/emc07r%dc%d_0",&row, &crystal);
896  }
897 
899  // Proto60
901  else if(node_path.Contains("Proto60")){
902  if(node_path.Contains("Passive") || node_path.Contains(TRegexp(".*PartAss_[0-9]*$")) ){
903  return kFALSE;
904  }
905  module=7;
906  copy=1;
907  Int_t type=1;
908 
909  if(node_path.Contains("CrystalType6aoPartAss"))
910  {
911  sscanf(node_path.Data(),"cave/Proto60_0/Active_1/Row%d_1/CrystalType6aoPartAss_%d/CrystalType6a_1",&row, &crystal);
912  } else if (node_path.Contains("CrystalType6boPartAss"))
913  {
914  sscanf(node_path.Data(),"cave/Proto60_0/Active_1/Row%d_1/CrystalType6boPartAss_%d/CrystalType6b_1",&row, &crystal);
915  type=2;
916  } else
917  {
918  cout<<"crystal name in Emc Proto: "<<node_path<<" missmatch pattern"<<endl;
919  return false;
920  }
921  crystal = ((crystal-1)%5)*2 +type;
922 
923 // printf("found node: %s\n Rowstring: %s crytalstring:%s\nrow: %d crystal: %d\n",node_path.Data(),((TObjString *)subStrL->At(1))->GetName(),((TObjString *)subStrL->At(2))->GetName(),row,crystal);
924  }
925 
927  //case of endcups of forward calorimeter from .dat file
929  else {
930  int tmp1;
931  sscanf(node_path.Data(),"cave/Emc%d_%d/emc%dr%dc%d_0",&tmp1, &copy, &module, &row, &crystal);
932  }
933  return true;
934 }
int row
Definition: anaLmdDigi.C:67
int col
Definition: anaLmdDigi.C:67
int getsign(const T &a)
Definition: PndEmc.cxx:53
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
static int next[96]
Definition: ranlxd.cxx:374
const mapper& PndEmcStructure::GetEmcX ( ) const
inline

Definition at line 35 of file PndEmcStructure.h.

References emcX.

Referenced by PndEmcApdHitProducer::Init(), and PndEmcHitProducer::Init().

35 { return emcX ;};
const mapper& PndEmcStructure::GetEmcY ( ) const
inline

Definition at line 36 of file PndEmcStructure.h.

References emcY.

Referenced by PndEmcApdHitProducer::Init(), and PndEmcHitProducer::Init().

36 { return emcY ;};
const mapper& PndEmcStructure::GetEmcZ ( ) const
inline

Definition at line 37 of file PndEmcStructure.h.

References emcZ.

Referenced by PndEmcApdHitProducer::Init(), and PndEmcHitProducer::Init().

37 { return emcZ ;};
const PndEmcTciXtalMap& PndEmcStructure::GetTciXtalMap ( ) const
inline
PndEmcStructure * PndEmcStructure::Instance ( )
static
PndEmcStructure * PndEmcStructure::Instance ( TGeoManager *  geoMan)
static

Definition at line 54 of file PndEmcStructure.cxx.

55 {
56  if (_instance == 0) {
57  if (geoMan==0){
58  cout<<"Empty geometry passed to PndEmcStructure"<<endl;
59  abort();
60  }
61  else
63  }
64  return _instance;
65 }
PndEmcStructure(TGeoManager *)
TGeoManager * geoMan
static PndEmcStructure * _instance
PndEmcTwoCoordIndex * PndEmcStructure::locateIndex ( double  theta,
double  phi 
) const

Definition at line 938 of file PndEmcStructure.cxx.

References fTciXtalMap, and vec.

Referenced by PndEmcClusterProperties::GravWhere(), PndEmcClusterProperties::LiloWhere(), PndEmcClusterProperties::LinearWhere(), and PndEmcClusterDistances::PndEmcClusterDistances().

939 {
940  // Find crystal with minimal angular difference with given direction
941  TVector3 vec(0,0,10);
942  vec.SetTheta(theta);
943  vec.SetPhi(phi);
944  PndEmcTwoCoordIndex *tci;
945  double diff=1000;
946 
947 // std::map <PndEmcTwoCoordIndex*, PndEmcXtal*>::const_iterator iter=fTciXtalMap.begin();
948  PndEmcTciXtalMap::const_iterator iter=fTciXtalMap.begin();
949  while(iter!=(fTciXtalMap).end())
950  {
951  TVector3 vec2= ((*iter).second)->frontCentre();
952  double tmpDiff=vec2.Angle(vec);
953 
954  if(tmpDiff<diff)
955  {
956  tci=(*iter).first;
957  diff=tmpDiff;
958  }
959 
960  iter++;
961  }
962 
963  return tci;
964 }
PndEmcTciXtalMap fTciXtalMap
stores crystal index coordinates (x,y) or (theta,phi)
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
void PndEmcStructure::Print ( string  filename,
Int_t  option = 1 
) const

Definition at line 968 of file PndEmcStructure.cxx.

References dPhi, dTheta, f, fTciXtalMap, PndEmcTwoCoordIndex::Index(), x, PndEmcTwoCoordIndex::XCoord(), y, PndEmcTwoCoordIndex::YCoord(), and z.

Referenced by MyMainFrame::DoRecreate(), and structure_analysis_emc().

969 {
970  ofstream f(filename.c_str());
971  //different option corresponds to different details level of output
972  if (option==1)
973  {
974  f<<"detID"<<"\t"<<"ThetaInd"<<"\t"<<"PhiInd"<<"\t"<<"theta_centre"<<"\t"<<"theta_frontface"<<"\t"<<"dTheta"<<"\t"<<"phi_centre"<<"\t"<<"phi_frontface"<<"\t"<<"dPhi"<<endl;
975  TVector3 centre, front_centre;
976  double theta_c, theta_f, phi_c, phi_f;
977  double dTheta, dPhi;
978 // std::map <PndEmcTwoCoordIndex*, PndEmcXtal*>::const_iterator iter=fTciXtalMap.begin();
979  PndEmcTciXtalMap::const_iterator iter=fTciXtalMap.begin();
980  PndEmcTwoCoordIndex *tci;
981  while(iter!=(fTciXtalMap).end())
982  {
983  tci = (*iter).first;
984  centre= ((*iter).second)->centre();
985  front_centre=((*iter).second)->frontCentre();
986  theta_c=centre.Theta()*TMath::RadToDeg();
987  phi_c=centre.Phi()*TMath::RadToDeg();
988  theta_f=front_centre.Theta()*TMath::RadToDeg();
989  phi_f=front_centre.Phi()*TMath::RadToDeg();
990 
991  dTheta=theta_c-theta_f;
992  dPhi=phi_c-phi_f;
993 
994  f<<tci->Index()<<"\t"<<tci->XCoord()<<"\t"<<tci->YCoord()<<"\t"<<theta_c<<"\t"<<theta_f<<"\t"<<dTheta<<"\t"<<phi_c<<"\t"<<phi_f<<"\t"<<dPhi<<endl;
995  iter++;
996  }
997  } else if (option==2)
998  {
999  f<<"detID"<<"\t"<<"ThetaInd"<<"\t"<<"PhiInd"<<"\t"<<"X"<<"\t"<<"Y"<<"\t"<<"Z"<<endl;
1000  TVector3 front_centre;
1001  TVector3 centre;
1002  double x, y, z;
1003 // std::map <PndEmcTwoCoordIndex*, PndEmcXtal*>::const_iterator iter=fTciXtalMap.begin();
1004  PndEmcTciXtalMap::const_iterator iter=fTciXtalMap.begin();
1005  PndEmcTwoCoordIndex *tci;
1006  while(iter!=(fTciXtalMap).end())
1007  {
1008  tci = (*iter).first;
1009  front_centre=((*iter).second)->frontCentre();
1010  centre=((*iter).second)->centre();
1011  //x=front_centre.X();y=front_centre.Y();z=front_centre.Z();
1012  x=centre.X();y=centre.Y();z=centre.Z();
1013 
1014  f<<tci->Index()<<"\t"<<tci->XCoord()<<"\t"<<tci->YCoord()<<"\t"<<x<<"\t"<<y<<"\t"<<z<<endl;
1015  iter++;
1016  }
1017 
1018  }
1019 
1020  f.close();
1021  cout<<"write emc structure to log file"<<endl;
1022 
1023 }
PndEmcTciXtalMap fTciXtalMap
stores crystal index coordinates (x,y) or (theta,phi)
TH1F * dTheta
Definition: anaLmdCluster.C:43
TFile * f
Definition: bump_analys.C:12
Double_t z
Double_t x
Double_t y
double dPhi
Definition: RiemannTest.C:17
const string filename

Member Data Documentation

PndEmcStructure * PndEmcStructure::_instance = 0
staticprivate

Definition at line 50 of file PndEmcStructure.h.

mapper PndEmcStructure::emcX
private

Definition at line 52 of file PndEmcStructure.h.

Referenced by GetEmcX(), and PndEmcStructure().

mapper PndEmcStructure::emcY
private

Definition at line 53 of file PndEmcStructure.h.

Referenced by GetEmcY(), and PndEmcStructure().

mapper PndEmcStructure::emcZ
private

Definition at line 54 of file PndEmcStructure.h.

Referenced by GetEmcZ(), and PndEmcStructure().

PndEmcTciXtalMap PndEmcStructure::fTciXtalMap
private

Definition at line 57 of file PndEmcStructure.h.

Referenced by GetTciXtalMap(), locateIndex(), PndEmcStructure(), and Print().


The documentation for this class was generated from the following files: