FairRoot/PandaRoot
Functions | Variables
createRootGeoFileBarrel_2018v1.C File Reference
#include "TFile.h"
#include "TGeoArb8.h"
#include "TGeoCompositeShape.h"
#include "TGeoManager.h"
#include "TGeoMaterial.h"
#include "TGeoMedium.h"
#include "TGeoMatrix.h"
#include "TGeoPcon.h"
#include "TGeoTube.h"
#include "TGeoVolume.h"
#include "TMath.h"
#include "TSystem.h"
#include "TVector3.h"
#include "TView.h"
#include "TVirtualPad.h"
#include "TTree.h"
#include <iostream>

Go to the source code of this file.

Functions

void geom ()
 
void ConstructSlice (TGeoVolume *vol)
 
void ConstructTargetSlice (TGeoVolume *vol)
 
void ConstructSuperModule (TGeoVolume *vol, Int_t id, Bool_t is_target=kFALSE)
 
void ConstructModule (TGeoVolume *vol, Int_t id, Int_t sign, Bool_t is_target=kFALSE)
 
void CalculateCrystalMatricesZ ()
 
void CalculateCrystalMatricesPhi ()
 
void CreateFrontInsertShapesAndMatricesZ ()
 
void CreateAlveoleShapesAndMatricesZ ()
 
TGeoTrap * create_trap (TString name, Int_t type, Double_t af, Double_t bf, Double_t cf, Double_t ar, Double_t br, Double_t cr, Double_t l)
 
TGeoTrap * create_trap (TString name, Int_t type, Double_t af, Double_t bf, Double_t cf, Double_t ar, Double_t br, Double_t cr, Double_t l, Int_t i)
 
void get_trap_vertices (TGeoTrap *, TVector3 *points)
 
void raytrace_x (TGeoManager *geom, TString filename)
 
void raytrace_z (TGeoManager *geom, TString filename)
 
int main ()
 

Variables

TGeoManager * geoMan
 
TGeoTrap * crystals_L [11]
 
TGeoTrap * crystals_R [11]
 
TGeoCompositeShape * wrappings_L [11]
 
TGeoCompositeShape * wrappings_R [11]
 
TGeoCompositeShape * front_inserts_L [11]
 
TGeoCompositeShape * front_inserts_R [11]
 
TGeoCompositeShape * alveoles_L_plus [43]
 
TGeoCompositeShape * alveoles_R_plus [43]
 
TGeoCompositeShape * alveoles_L_minus [28]
 
TGeoCompositeShape * alveoles_R_minus [28]
 
TGeoCombiTrans * matrices_minus_left [28]
 
TGeoCombiTrans * matrices_minus_right [28]
 
TGeoCombiTrans * matrices_plus_left [43]
 
TGeoCombiTrans * matrices_plus_right [43]
 
TGeoTranslation * matrices_fe2cry_left [11]
 
TGeoTranslation * matrices_fe2cry_right [11]
 
TGeoMatrix * matrices_alveole_minus_left [28]
 
TGeoMatrix * matrices_alveole_minus_right [28]
 
TGeoMatrix * matrices_alveole_plus_left [43]
 
TGeoMatrix * matrices_alveole_plus_right [43]
 
Double_t crystal_dz [43]
 
Double_t crystal_theta [44]
 
TGeoMatrix * matrices_phi [5]
 
TGeoTranslation * matrices_crystal2alveole_L [11]
 
TGeoTranslation * matrices_crystal2alveole_R [11]
 
Double_t AF [11] = {2.121, 2.118, 2.117, 2.117, 2.117, 2.119, 2.122, 2.123, 2.123, 2.125, 2.125}
 
Double_t BF [11] = {2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128}
 
Double_t CF [11] = {2.127, 2.139, 2.151, 2.160, 2.169, 2.178, 2.186, 2.191, 2.195, 2.200, 2.202}
 
Double_t AR [11] = {2.904, 2.878, 2.836, 2.790, 2.735, 2.672, 2.623, 2.570, 2.514, 2.470, 2.435}
 
Double_t BR [11] = {2.875, 2.875, 2.875, 2.722, 2.722, 2.722, 2.547, 2.547, 2.547, 2.442, 2.442}
 
Double_t CR [11] = {2.912, 2.907, 2.881, 2.845, 2.801, 2.747, 2.699, 2.651, 2.600, 2.556, 2.523}
 
Double_t theta_hole_corr = -(4.997 + 3.862213309 + 90.)
 

Function Documentation

void CalculateCrystalMatricesPhi ( )

Definition at line 681 of file createRootGeoFileBarrel_2018v1.C.

References Double_t, i, matrices_phi, rot, rot1, theta, and theta_hole_corr.

Referenced by geom().

681  {
682  Double_t alpha_alveole[] = {3.862213309, 3.945818621, 4.028979744, 4.11168967, 4.193941577};
683  Double_t theta_alveole[] = {0, 4.416394688, 4.416838877, 4.417290073, 4.417748093};
685  for (Int_t i = 0; i < 5; i++) {
686  theta += theta_alveole[i];
687  TGeoRotation rot(TString::Format("rot_alveole_row_%d", i),
688  90, 360-theta, 90, 90-theta, 0, 0); // clockwise
689 
690  TGeoTranslation tra1(-57., 0, 0);
691  TGeoRotation rot1("rot1_row2", 90, 360-alpha_alveole[i], 90, 90-alpha_alveole[i], 0, 0);
692  TGeoTranslation tra2(57., 0, 0);
693  matrices_phi[i] = new TGeoHMatrix(rot * tra2 * rot1 * tra1);
694  matrices_phi[i]->SetName(TString::Format("mat_phi_%d", i));
695  }
696 }
Int_t i
Definition: run_full.C:25
TGeoRotation * rot1
Double_t theta_hole_corr
Double_t
TGeoMatrix * matrices_phi[5]
TGeoRotation rot
void CalculateCrystalMatricesZ ( )

Definition at line 507 of file createRootGeoFileBarrel_2018v1.C.

References CAMath::Cos(), crystal_dz, crystal_theta, crystals_L, crystals_R, Double_t, dz, get_trap_vertices(), h1, h2, i, matrices_minus_left, matrices_minus_right, matrices_plus_left, matrices_plus_right, Pi, and theta.

Referenced by geom().

507  {
508  for (Int_t type = 0; type < 4; type++) { // type: 0 left-, 1 right+, 2 right-, 3 left+
509  const Double_t crystal_gap = 0.068;
510  const Double_t alveole_gap = 0.09;
511  Double_t module_gap[] = {0.12, 0.25, 0.3, 0.33}; // half of the first gap
512 
513  Double_t trans_z = 0;
514  Double_t trans_theta = 0;
515 
516  Int_t nModule;
517  if (type == 0 || type == 2) nModule = 3;
518  else nModule = 4;
519  Int_t i = 0;
520  crystal_theta[0] = 0.;
521  for (Int_t iModule = 0; iModule < nModule; iModule++) {
522  for (Int_t iAlveole = 0; iAlveole < 3; iAlveole++) {
523  for (Int_t iCrystal = 0; iCrystal < 4; iCrystal++) {
524  Int_t module_id;
525  if (type == 0 || type == 2) module_id = 3 - iModule;
526  else module_id = iModule + 4;
527  Int_t alveole_id = iModule * 3 + iAlveole + 1;
528  Int_t crystal_id = iCrystal + 1;
529 
530  if (module_id == 1 && iAlveole > 0) break; // SuperModule1 only have 1 alveole
531  if (module_id == 7 && iAlveole > 1) break; // Module7 only have 2 alveole
532  if (module_id == 7 && iAlveole == 1 && iCrystal > 2) break; // The last alveole has 3 crystal colomns
533 
534  // calculate crystal type
535  TGeoTrap* crystal;
536  if (type == 0 || type == 3) crystal = crystals_L[alveole_id-1];
537  else crystal = crystals_R[alveole_id-1];
538 
539  Double_t gap = crystal_gap;
540  if (iCrystal == 0) {
541  gap = alveole_gap;
542  if (iAlveole == 0) gap = module_gap[iModule];
543  }
544 
545  Double_t bf = crystal->GetH1()*2;
546  Double_t br = crystal->GetH2()*2;
547  Double_t l = crystal->GetDz()*2;
548  Double_t theta = TMath::ATan((br-bf)/l);
549 
550  Double_t bl1 = crystal->GetBl1();
551  Double_t bl2 = crystal->GetBl2();
552  Double_t tl1 = crystal->GetTl1();
553  Double_t tl2 = crystal->GetTl2();
554  Double_t h1 = crystal->GetH1();
555  Double_t h2 = crystal->GetH2();
556  Double_t dz = crystal->GetDz();
557 
558  // calculate the center of the crystal
559  TVector3 center2pivot;
560  TVector3 vertices[8];
561  get_trap_vertices(crystal, vertices);
562  if (type == 0) {
563  trans_z -= (gap + bf) / TMath::Cos(trans_theta); // trans_z < 0
564  center2pivot = -vertices[1]; // rotate c2p from ROOT coordinate to experimental coordinate
565  }
566  else if (type == 1) {
567  trans_z += (gap + bf) / TMath::Cos(trans_theta); // trans_z < 0
568  center2pivot = -vertices[2]; // rotate c2p from ROOT coordinate to experimental coordinate
569  }
570  else if (type == 2) {
571  trans_z -= (gap + bf) / TMath::Cos(trans_theta); // trans_z < 0
572  center2pivot = -vertices[2];
573  }
574  else {
575  trans_z += (gap + bf) / TMath::Cos(trans_theta); // trans_z < 0
576  center2pivot = -vertices[1];
577  }
578  if (type == 0) {
579  center2pivot.RotateX(-TMath::Pi()/2.);
580  center2pivot.RotateZ(-TMath::Pi()/2.);
581  }
582  else if (type == 1) {
583  center2pivot.RotateZ(TMath::Pi());
584  center2pivot.RotateX(-TMath::Pi()/2.);
585  center2pivot.RotateZ(-TMath::Pi()/2.);
586  }
587  else if (type == 2) {
588  center2pivot.RotateZ(TMath::Pi());
589  center2pivot.RotateX(TMath::Pi()/2.);
590  center2pivot.RotateZ(TMath::Pi()/2.);
591  }
592  else {
593  center2pivot.RotateX(TMath::Pi()/2.);
594  center2pivot.RotateZ(TMath::Pi()/2.);
595  }
596  center2pivot.RotateY(trans_theta);
597  TVector3 pivot;
598  if (type == 0 || type == 1)
599  pivot = TVector3(57.+0.112126393, -0.068, trans_z + 3.7); // EMC radius is 570mm, crystal offset w.r.t. target is 37mm
600  else
601  pivot = TVector3(57., 0, trans_z + 3.7);
602  TVector3 center = pivot + center2pivot;
603 
604  // calculate the rotation matrix and place the crystal
605  Double_t trans_theta_deg = trans_theta/TMath::Pi()*180.;
606  TGeoRotation* crystal_rot;
607  char avleole_sign;
608  char crystal_type;
609  if (type == 0) {
610  avleole_sign = 'm';
611  crystal_type = 'L';
612  }
613  else if (type == 1) {
614  avleole_sign = 'p';
615  crystal_type = 'R';
616  }
617  else if (type == 1) {
618  avleole_sign = 'm';
619  crystal_type = 'R';
620  }
621  else {
622  avleole_sign = 'p';
623  crystal_type = 'L';
624  }
625  if (type == 0) {
626  crystal_rot = new TGeoRotation(
627  TString::Format("rot-%d%c-%c%d",
628  alveole_id, avleole_sign, crystal_type, crystal_id),
629  90, 270, 180-trans_theta_deg, 180, 90+trans_theta_deg, 0);
630  }
631  else if (type == 1) {
632  crystal_rot = new TGeoRotation(
633  TString::Format("rot-%d%c-%c%d",
634  alveole_id, avleole_sign, crystal_type, crystal_id),
635  90, 90, trans_theta_deg, 0, 90+trans_theta_deg, 0);
636  }
637  else if (type == 2) {
638  crystal_rot = new TGeoRotation(
639  TString::Format("rot-%d%c-%c%d",
640  alveole_id, avleole_sign, crystal_type, crystal_id),
641  90, 270, 180-trans_theta_deg, 180, 90+trans_theta_deg, 0);
642  }
643  else {
644  crystal_rot = new TGeoRotation(
645  TString::Format("rot-%d%c-%c%d",
646  alveole_id, avleole_sign, crystal_type, crystal_id),
647  90, 90, trans_theta_deg, 0, 90+trans_theta_deg, 0);
648  }
649  crystal_rot->RegisterYourself();
650 
651  TGeoCombiTrans* crystal_ctrans = new TGeoCombiTrans(TString::Format("trans-%d%c-%c%d",
652  alveole_id, avleole_sign, crystal_type, crystal_id),
653  center.X(), center.Y(), center.Z(), crystal_rot);
654 
655  if (type == 0)
656  matrices_minus_left[i] = crystal_ctrans;
657  else if (type == 1)
658  matrices_plus_right[i] = crystal_ctrans;
659  else if (type == 2)
660  matrices_minus_right[i] = crystal_ctrans;
661  else
662  matrices_plus_left[i] = crystal_ctrans;
663 
664  if (type == 0 || type == 2)
665  trans_theta += theta; // trans_theta > 0
666  else
667  trans_theta -= theta;
668 
669  if (type == 1) {
670  crystal_theta[i+1] = -trans_theta;
671  crystal_dz[i] = trans_z;
672  }
673 
674  i++;
675  } // crystal
676  } // module (alveole)
677  } // supermodule
678  } // type
679 }
void get_trap_vertices(TGeoTrap *, TVector3 *points)
TGeoCombiTrans * matrices_plus_left[43]
TGeoTrap * crystals_L[11]
Int_t i
Definition: run_full.C:25
static T Cos(const T &x)
Definition: PndCAMath.h:43
Double_t
TGeoCombiTrans * matrices_minus_right[28]
TGeoCombiTrans * matrices_plus_right[43]
Double_t crystal_dz[43]
TGeoCombiTrans * matrices_minus_left[28]
Double_t crystal_theta[44]
Double_t Pi
TGeoTrap * crystals_R[11]
void ConstructModule ( TGeoVolume *  vol,
Int_t  id,
Int_t  sign,
Bool_t  is_target = kFALSE 
)

Definition at line 386 of file createRootGeoFileBarrel_2018v1.C.

References alveoles_L_minus, alveoles_L_plus, alveoles_R_minus, alveoles_R_plus, crystals_L, crystals_R, front_inserts_L, front_inserts_R, geoMan, i, matrices_alveole_minus_left, matrices_alveole_minus_right, matrices_alveole_plus_left, matrices_alveole_plus_right, matrices_fe2cry_left, matrices_fe2cry_right, matrices_minus_left, matrices_minus_right, matrices_phi, matrices_plus_left, matrices_plus_right, TString, wrappings_L, and wrappings_R.

Referenced by ConstructSuperModule().

386  {
387  // Volumes
388  TGeoVolume* Crystal1;
389  TGeoVolume* Crystal2;
390  TGeoVolume* Wrapping1;
391  TGeoVolume* Wrapping2;
392  TGeoVolume* FrontInsert1;
393  TGeoVolume* FrontInsert2;
394  TGeoVolume* Alveole1;
395  TGeoVolume* Alveole2;
396  //TGeoMaterial *matVacuum = new TGeoMaterial("Vacuum", 0, 0, 0);
397  //TGeoMedium *Vacuum = new TGeoMedium("Vacuum",0, matVacuum);
398 
399  Int_t nCrystals = (id == 11) ? 3 : 4;
400  char char_sign = (sign > 0) ? 'p' : 'm';
401  char char_type1 = (sign > 0) ? 'R' : 'L';
402  char char_type2 = (sign > 0) ? 'L' : 'R';
403 
404  Int_t j_start = (id-1) * 4;
405  Int_t j_end = j_start + nCrystals;
406  for (Int_t i = 0; i < 5; i++) {
407  for(Int_t j = j_start; j < j_end; j++) {
408  if (is_target) {
409  if (id == 1 && (0 < i && i < 4)) continue;
410  }
411 
412  TString crystal_name1 = TString::Format("Crystal-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1);
413  TString crystal_name2 = TString::Format("Crystal-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1);
414  TString wrapping_name1 = TString::Format("Wrapping-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1);
415  TString wrapping_name2 = TString::Format("Wrapping-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1);
416  TString front_insert_name1 = TString::Format("FrontInsert-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1);
417  TString front_insert_name2 = TString::Format("FrontInsert-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1);
418  TString alveole_name1 = TString::Format("Alveole-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1);
419  TString alveole_name2 = TString::Format("Alveole-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1);
420  if (sign < 0) {
421  Crystal1 = new TGeoVolume(crystal_name1, crystals_L[id-1], geoMan->GetMedium("PWO"));
422  Crystal2 = new TGeoVolume(crystal_name2, crystals_R[id-1], geoMan->GetMedium("PWO"));
423  Wrapping1 = new TGeoVolume(wrapping_name1, wrappings_L[id-1], geoMan->GetMedium("VM2000"));
424  Wrapping2 = new TGeoVolume(wrapping_name2, wrappings_R[id-1], geoMan->GetMedium("VM2000"));
425  FrontInsert1 = new TGeoVolume(front_insert_name1, front_inserts_L[id-1], geoMan->GetMedium("ABS"));
426  FrontInsert2 = new TGeoVolume(front_insert_name2, front_inserts_R[id-1], geoMan->GetMedium("ABS"));
427  Alveole1 = new TGeoVolume(alveole_name1, alveoles_L_minus[j], geoMan->GetMedium("Prepreg"));
428  Alveole2 = new TGeoVolume(alveole_name2, alveoles_R_minus[j], geoMan->GetMedium("Prepreg"));
429  }
430  else {
431  Crystal1 = new TGeoVolume(crystal_name1, crystals_R[id-1], geoMan->GetMedium("PWO"));
432  Crystal2 = new TGeoVolume(crystal_name2, crystals_L[id-1], geoMan->GetMedium("PWO"));
433  Wrapping1 = new TGeoVolume(wrapping_name1, wrappings_R[id-1], geoMan->GetMedium("VM2000"));
434  Wrapping2 = new TGeoVolume(wrapping_name2, wrappings_L[id-1], geoMan->GetMedium("VM2000"));
435  FrontInsert1 = new TGeoVolume(front_insert_name1, front_inserts_R[id-1], geoMan->GetMedium("ABS"));
436  FrontInsert2 = new TGeoVolume(front_insert_name2, front_inserts_L[id-1], geoMan->GetMedium("ABS"));
437  Alveole1 = new TGeoVolume(alveole_name1, alveoles_R_plus[j], geoMan->GetMedium("Prepreg"));
438  Alveole2 = new TGeoVolume(alveole_name2, alveoles_L_plus[j], geoMan->GetMedium("Prepreg"));
439  }
440  Crystal1->SetLineColor(kCyan);
441  Crystal2->SetLineColor(kCyan);
442  Wrapping1->SetLineColor(kViolet);
443  Wrapping2->SetLineColor(kViolet);
444  FrontInsert1->SetLineColor(kWhite);
445  FrontInsert2->SetLineColor(kWhite);
446  Alveole1->SetLineColor(kYellow);
447  Alveole2->SetLineColor(kYellow);
448 
449  // Matrices
450  TGeoHMatrix* mat1_crystal;
451  TGeoHMatrix* mat2_crystal;
452  if (sign < 0) {
453  mat1_crystal = new TGeoHMatrix((*matrices_phi[i])*(*matrices_minus_left[j]));
454  mat2_crystal = new TGeoHMatrix((*matrices_phi[i])*(*matrices_minus_right[j]));
455  }
456  else {
457  mat1_crystal = new TGeoHMatrix((*matrices_phi[i])*(*matrices_plus_right[j]));
458  mat2_crystal = new TGeoHMatrix((*matrices_phi[i])*(*matrices_plus_left[j]));
459  }
460  mat1_crystal->SetName(TString::Format("mat-crystal-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1));
461  mat2_crystal->SetName(TString::Format("mat-crystal-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1));
462 
463  TGeoHMatrix* mat1_front_insert;
464  TGeoHMatrix* mat2_front_insert;
465  if (sign < 0) {
466  mat1_front_insert = new TGeoHMatrix((*matrices_phi[i])*(*matrices_minus_left[j])*(*matrices_fe2cry_left[id-1]));
467  mat2_front_insert = new TGeoHMatrix((*matrices_phi[i])*(*matrices_minus_right[j])*(*matrices_fe2cry_right[id-1]));
468  }
469  else {
470  mat1_front_insert = new TGeoHMatrix((*matrices_phi[i])*(*matrices_plus_right[j])*(*matrices_fe2cry_right[id-1]));
471  mat2_front_insert = new TGeoHMatrix((*matrices_phi[i])*(*matrices_plus_left[j])*(*matrices_fe2cry_left[id-1]));
472  }
473  mat1_front_insert->SetName(TString::Format("mat-front-insert-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1));
474  mat2_front_insert->SetName(TString::Format("mat-front-insert-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1));
475 
476  TGeoHMatrix* mat1_alveole;
477  TGeoHMatrix* mat2_alveole;
478  if (sign < 0) {
479  mat1_alveole = new TGeoHMatrix((*matrices_phi[i])*
481  mat2_alveole = new TGeoHMatrix((*matrices_phi[i])*
483  }
484  else {
485  mat1_alveole = new TGeoHMatrix((*matrices_phi[i])*
487  mat2_alveole = new TGeoHMatrix((*matrices_phi[i])*
489  }
490  mat1_alveole->SetName(TString::Format("mat-alveole-%d%c-%c%d", id, char_sign, char_type1, j-j_start+1));
491  mat2_alveole->SetName(TString::Format("mat-alveole-%d%c-%c%d", id, char_sign, char_type2, j-j_start+1));
492 
493 
494  // Nodes
495  vol->AddNode(Alveole1, i+1, mat1_alveole);
496  vol->AddNode(Alveole2, i+1, mat2_alveole);
497  vol->AddNode(FrontInsert1, i+1, mat1_front_insert);
498  vol->AddNode(FrontInsert2, i+1, mat2_front_insert);
499  vol->AddNode(Wrapping1, i+1, mat1_crystal);
500  vol->AddNode(Wrapping2, i+1, mat2_crystal);
501  vol->AddNode(Crystal1, i+1, mat1_crystal);
502  vol->AddNode(Crystal2, i+1, mat2_crystal);
503  }
504  }
505 }
TGeoMatrix * matrices_alveole_minus_right[28]
TGeoTranslation * matrices_fe2cry_left[11]
TGeoCombiTrans * matrices_plus_left[43]
TGeoMatrix * matrices_alveole_plus_right[43]
TGeoTrap * crystals_L[11]
Int_t i
Definition: run_full.C:25
TGeoMatrix * matrices_alveole_plus_left[43]
TGeoManager * geoMan
TGeoCompositeShape * wrappings_R[11]
TGeoMatrix * matrices_alveole_minus_left[28]
TGeoCompositeShape * alveoles_R_plus[43]
TGeoCompositeShape * wrappings_L[11]
TGeoCompositeShape * alveoles_L_minus[28]
TGeoCompositeShape * front_inserts_L[11]
TGeoCompositeShape * alveoles_R_minus[28]
TGeoCombiTrans * matrices_minus_right[28]
TGeoCompositeShape * front_inserts_R[11]
TGeoMatrix * matrices_phi[5]
TGeoTranslation * matrices_fe2cry_right[11]
TGeoCombiTrans * matrices_plus_right[43]
TGeoCompositeShape * alveoles_L_plus[43]
int sign(T val)
Definition: PndCADef.h:48
TGeoCombiTrans * matrices_minus_left[28]
TGeoTrap * crystals_R[11]
void ConstructSlice ( TGeoVolume *  vol)

Definition at line 210 of file createRootGeoFileBarrel_2018v1.C.

References ConstructSuperModule(), geoMan, and theta_hole_corr.

Referenced by geom().

210  {
211  TGeoVolumeAssembly* SuperModule1 = new TGeoVolumeAssembly("SuperModule1");
212  TGeoVolumeAssembly* SuperModule2 = new TGeoVolumeAssembly("SuperModule2");
213  TGeoVolumeAssembly* SuperModule3 = new TGeoVolumeAssembly("SuperModule3");
214  TGeoVolumeAssembly* SuperModule4 = new TGeoVolumeAssembly("SuperModule4");
215  TGeoVolumeAssembly* SuperModule5 = new TGeoVolumeAssembly("SuperModule5");
216  TGeoVolumeAssembly* SuperModule6 = new TGeoVolumeAssembly("SuperModule6");
217  TGeoVolumeAssembly* SuperModule7 = new TGeoVolumeAssembly("SuperModule7");
218 
219  ConstructSuperModule(SuperModule1, 1);
220  ConstructSuperModule(SuperModule2, 2);
221  ConstructSuperModule(SuperModule3, 3);
222  ConstructSuperModule(SuperModule4, 4);
223  ConstructSuperModule(SuperModule5, 5);
224  ConstructSuperModule(SuperModule6, 6);
225  ConstructSuperModule(SuperModule7, 7);
226 
227  vol->AddNode(SuperModule1, 1);
228  vol->AddNode(SuperModule2, 1);
229  vol->AddNode(SuperModule3, 1);
230  vol->AddNode(SuperModule4, 1);
231  vol->AddNode(SuperModule5, 1);
232  vol->AddNode(SuperModule6, 1);
233  vol->AddNode(SuperModule7, 1);
234 
235  // Simplified front end structures
236  TGeoTubeSeg* front_logement = new TGeoTubeSeg(55.0, 55.1, 106.6, -theta_hole_corr-19.89, -theta_hole_corr+2.21);
237  TGeoVolume* FrontLogement = new TGeoVolume("Front-Logement", front_logement, geoMan->GetMedium("Polyurethane"));
238  FrontLogement->SetLineColor(kBlue);
239  TGeoTranslation* mat_logement = new TGeoTranslation("mat-Front-Logement", 0., 0., 36.9);
240  vol->AddNode(FrontLogement, 1, mat_logement);
241 
242  TGeoTubeSeg* front_al_plate = new TGeoTubeSeg(54.2, 54.3, 106.6, -theta_hole_corr-19.89, -theta_hole_corr+2.21);
243  TGeoVolume* FrontAlPlate = new TGeoVolume("Front-Al-Plate", front_al_plate, geoMan->GetMedium("Aluminum"));
244  FrontAlPlate->SetLineColor(kDarkTerrain);
245  TGeoTranslation* mat_al_plate = new TGeoTranslation("mat-Front-Al-Plate", 0., 0., 36.9);
246  vol->AddNode(FrontAlPlate, 1, mat_al_plate);
247 }
TGeoManager * geoMan
Double_t theta_hole_corr
void ConstructSuperModule(TGeoVolume *vol, Int_t id, Bool_t is_target=kFALSE)
void ConstructSuperModule ( TGeoVolume *  vol,
Int_t  id,
Bool_t  is_target = kFALSE 
)

Definition at line 304 of file createRootGeoFileBarrel_2018v1.C.

References ConstructModule().

Referenced by ConstructSlice(), and ConstructTargetSlice().

304  {
305  if (id == 1) {
306  TGeoVolumeAssembly* Module7m = new TGeoVolumeAssembly("Module7-");
307  ConstructModule(Module7m, 7, -1, is_target);
308  vol->AddNode(Module7m, 1);
309  }
310  else if (id == 2) {
311  TGeoVolumeAssembly* Module6m = new TGeoVolumeAssembly("Module6-");
312  ConstructModule(Module6m, 6, -1, is_target);
313  vol->AddNode(Module6m, 1);
314 
315  TGeoVolumeAssembly* Module5m = new TGeoVolumeAssembly("Module5-");
316  ConstructModule(Module5m, 5, -1, is_target);
317  vol->AddNode(Module5m, 1);
318 
319  TGeoVolumeAssembly* Module4m = new TGeoVolumeAssembly("Module4-");
320  ConstructModule(Module4m, 4, -1, is_target);
321  vol->AddNode(Module4m, 1);
322  }
323  else if (id == 3) {
324  TGeoVolumeAssembly* Module3m = new TGeoVolumeAssembly("Module3-");
325  ConstructModule(Module3m, 3, -1, is_target);
326  vol->AddNode(Module3m, 1);
327 
328  TGeoVolumeAssembly* Module2m = new TGeoVolumeAssembly("Module2-");
329  ConstructModule(Module2m, 2, -1, is_target);
330  vol->AddNode(Module2m, 1);
331 
332  TGeoVolumeAssembly* Module1m = new TGeoVolumeAssembly("Module1-");
333  ConstructModule(Module1m, 1, -1, is_target);
334  vol->AddNode(Module1m, 1);
335  }
336  else if (id == 4) {
337  TGeoVolumeAssembly* Module1p = new TGeoVolumeAssembly("Module1+");
338  ConstructModule(Module1p, 1, +1, is_target);
339  vol->AddNode(Module1p, 1);
340 
341  TGeoVolumeAssembly* Module2p = new TGeoVolumeAssembly("Module2+");
342  ConstructModule(Module2p, 2, +1, is_target);
343  vol->AddNode(Module2p, 1);
344 
345  TGeoVolumeAssembly* Module3p = new TGeoVolumeAssembly("Module3+");
346  ConstructModule(Module3p, 3, +1, is_target);
347  vol->AddNode(Module3p, 1);
348  }
349  else if (id == 5) {
350  TGeoVolumeAssembly* Module4p = new TGeoVolumeAssembly("Module4+");
351  ConstructModule(Module4p, 4, +1, is_target);
352  vol->AddNode(Module4p, 1);
353 
354  TGeoVolumeAssembly* Module5p = new TGeoVolumeAssembly("Module5+");
355  ConstructModule(Module5p, 5, +1, is_target);
356  vol->AddNode(Module5p, 1);
357 
358  TGeoVolumeAssembly* Module6p = new TGeoVolumeAssembly("Module6+");
359  ConstructModule(Module6p, 6, +1, is_target);
360  vol->AddNode(Module6p, 1);
361  }
362  else if (id == 6) {
363  TGeoVolumeAssembly* Module7p = new TGeoVolumeAssembly("Module7+");
364  ConstructModule(Module7p, 7, +1, is_target);
365  vol->AddNode(Module7p, 1);
366 
367  TGeoVolumeAssembly* Module8p = new TGeoVolumeAssembly("Module8+");
368  ConstructModule(Module8p, 8, +1, is_target);
369  vol->AddNode(Module8p, 1);
370 
371  TGeoVolumeAssembly* Module9p = new TGeoVolumeAssembly("Module9+");
372  ConstructModule(Module9p, 9, +1, is_target);
373  vol->AddNode(Module9p, 1);
374  }
375  else if (id == 7) {
376  TGeoVolumeAssembly* Module10p = new TGeoVolumeAssembly("Module10+");
377  ConstructModule(Module10p, 10, +1, is_target);
378  vol->AddNode(Module10p, 1);
379 
380  TGeoVolumeAssembly* Module11p = new TGeoVolumeAssembly("Module11+");
381  ConstructModule(Module11p, 11, +1, is_target);
382  vol->AddNode(Module11p, 1);
383  }
384 }
void ConstructModule(TGeoVolume *vol, Int_t id, Int_t sign, Bool_t is_target=kFALSE)
void ConstructTargetSlice ( TGeoVolume *  vol)

Definition at line 249 of file createRootGeoFileBarrel_2018v1.C.

References ConstructSuperModule(), geoMan, and theta_hole_corr.

Referenced by geom().

249  {
250  TGeoVolumeAssembly* SuperModule1 = new TGeoVolumeAssembly("SuperModule1_Target");
251  TGeoVolumeAssembly* SuperModule2 = new TGeoVolumeAssembly("SuperModule2_Target");
252  TGeoVolumeAssembly* SuperModule3 = new TGeoVolumeAssembly("SuperModule3_Target");
253  TGeoVolumeAssembly* SuperModule4 = new TGeoVolumeAssembly("SuperModule4_Target");
254  TGeoVolumeAssembly* SuperModule5 = new TGeoVolumeAssembly("SuperModule5_Target");
255  TGeoVolumeAssembly* SuperModule6 = new TGeoVolumeAssembly("SuperModule6_Target");
256  TGeoVolumeAssembly* SuperModule7 = new TGeoVolumeAssembly("SuperModule7_Target");
257 
258  ConstructSuperModule(SuperModule1, 1, kTRUE);
259  ConstructSuperModule(SuperModule2, 2, kTRUE);
260  ConstructSuperModule(SuperModule3, 3, kTRUE);
261  ConstructSuperModule(SuperModule4, 4, kTRUE);
262  ConstructSuperModule(SuperModule5, 5, kTRUE);
263  ConstructSuperModule(SuperModule6, 6, kTRUE);
264  ConstructSuperModule(SuperModule7, 7, kTRUE);
265 
266  vol->AddNode(SuperModule1, 1);
267  vol->AddNode(SuperModule2, 1);
268  vol->AddNode(SuperModule3, 1);
269  vol->AddNode(SuperModule4, 1);
270  vol->AddNode(SuperModule5, 1);
271  vol->AddNode(SuperModule6, 1);
272  vol->AddNode(SuperModule7, 1);
273 
274  /* Simplified front end structures */
275  // Logement for cooling tubes
276  TGeoTubeSeg* front_logement = new TGeoTubeSeg("front_logement", 55.0, 55.1, 106.6, -theta_hole_corr-19.89, -theta_hole_corr+2.21);
277  TGeoTranslation* traLogement = new TGeoTranslation("traLogement", 0., 0., 36.9);
278  traLogement->RegisterYourself();
279  TGeoTube* holeLogement = new TGeoTube("holeLogement", 0., 7.6, 100.);
280  TGeoRotation* rotHoleLogement = new TGeoRotation("rotHoleLogement", 90, 0, 0, 0, 90, 270);
281  TGeoCombiTrans* matHoleLogement = new TGeoCombiTrans("matHoleLogement", 0., 0., 3.7, rotHoleLogement);
282  matHoleLogement->RegisterYourself();
283  TGeoCompositeShape* front_logement_target =
284  new TGeoCompositeShape("front_logement_target", "front_logement:traLogement - holeLogement:matHoleLogement");
285  TGeoVolume* FrontLogement = new TGeoVolume("Front-Logement-Target", front_logement_target, geoMan->GetMedium("Polyurethane"));
286  FrontLogement->SetLineColor(kBlue);
287  vol->AddNode(FrontLogement, 1);
288 
289  // Aluminum plate for vacuum
290  TGeoTubeSeg* front_al_plate = new TGeoTubeSeg("front_al_plate", 54.2, 54.3, 106.6, -theta_hole_corr-19.89, -theta_hole_corr+2.21);
291  TGeoTranslation* traAlPlate = new TGeoTranslation("traAlPlate", 0., 0., 36.9);
292  traAlPlate->RegisterYourself();
293  TGeoTube* holeAlPlate = new TGeoTube("holeAlPlate", 0., 7.6, 100.);
294  TGeoRotation* rotHoleAlPlate = new TGeoRotation("rotHoleAlPlate", 90, 0, 0, 0, 90, 270);
295  TGeoCombiTrans* matHoleAlPlate = new TGeoCombiTrans("matHoleAlPlate", 0., 0., 3.7, rotHoleAlPlate);
296  matHoleAlPlate->RegisterYourself();
297  TGeoCompositeShape* front_al_plate_target =
298  new TGeoCompositeShape("front_al_plate_target", "front_al_plate:traAlPlate - holeAlPlate:matHoleAlPlate");
299  TGeoVolume* FrontAlPlate = new TGeoVolume("Front-Al-Plate-Target", front_al_plate_target, geoMan->GetMedium("Aluminum"));
300  FrontAlPlate->SetLineColor(kDarkTerrain);
301  vol->AddNode(FrontAlPlate, 1);
302 }
TGeoManager * geoMan
Double_t theta_hole_corr
void ConstructSuperModule(TGeoVolume *vol, Int_t id, Bool_t is_target=kFALSE)
TGeoTrap * create_trap ( TString  name,
Int_t  type,
Double_t  af,
Double_t  bf,
Double_t  cf,
Double_t  ar,
Double_t  br,
Double_t  cr,
Double_t  l 
)

Definition at line 1110 of file createRootGeoFileBarrel_2018v1.C.

References Double_t, dz, h1, h2, phi, Pi, and theta.

Referenced by CreateAlveoleShapesAndMatricesZ(), CreateFrontInsertShapesAndMatricesZ(), and geom().

1112 {
1113  Double_t x1 = (af + cf) / 4;
1114  Double_t y1 = bf / 2;
1115  Double_t z1 = 0.;
1116  Double_t x2 = (ar + cr) / 4;
1117  Double_t y2 = br / 2;
1118  Double_t z2 = l;
1119  Double_t x21 = x2 - x1;
1120  Double_t y21 = y2 - y1;
1121  Double_t z21 = z2 - z1;
1122  TVector3 v21(x21, y21, z21);
1123 
1124  Double_t dz = l / 2.;
1125  Double_t theta = v21.Theta() / TMath::Pi() * 180.;
1126  Double_t phi = v21.Phi() / TMath::Pi() * 180.;
1127  Double_t h1 = bf / 2;
1128  Double_t bl1 = cf / 2;
1129  Double_t tl1 = af / 2;
1130  Double_t alpha1 = -TMath::ATan((cf-af)/2/bf) / TMath::Pi() * 180.;
1131  Double_t h2 = br / 2;
1132  Double_t bl2 = cr / 2;
1133  Double_t tl2 = ar / 2;
1134  Double_t alpha2 = -TMath::ATan((cr-ar)/2/br) / TMath::Pi() * 180.;
1135 
1136  TGeoTrap* trap;
1137  if (type > 0)
1138  trap = new TGeoTrap(name, dz, theta, phi, h1, bl1, tl1, alpha1,
1139  h2, bl2, tl2, alpha2);
1140  else
1141  trap = new TGeoTrap(name, dz, theta, 180-phi, h1, bl1, tl1, -alpha1,
1142  h2, bl2, tl2, -alpha2);
1143 
1144  return trap;
1145 }
Double_t
TString name
Double_t Pi
TGeoTrap * create_trap ( TString  name,
Int_t  type,
Double_t  af,
Double_t  bf,
Double_t  cf,
Double_t  ar,
Double_t  br,
Double_t  cr,
Double_t  l,
Int_t  i 
)

Definition at line 1147 of file createRootGeoFileBarrel_2018v1.C.

References crystal_theta, Double_t, dz, h1, h2, phi, Pi, CAMath::Tan(), and theta.

1149 {
1150  Double_t x1 = (af + cf) / 4;
1151  Double_t y1 = bf / 2;
1152  Double_t z1 = 0.;
1153  Double_t x2 = (ar + cr) / 4;
1154  Double_t y2 = br / 2 + l * TMath::Tan(crystal_theta[i]);
1155  Double_t z2 = l;
1156  Double_t x21 = x2 - x1;
1157  Double_t y21 = y2 - y1;
1158  Double_t z21 = z2 - z1;
1159  TVector3 v21(x21, y21, z21);
1160 
1161  Double_t dz = l / 2.;
1162  Double_t theta = v21.Theta() / TMath::Pi() * 180.;
1163  Double_t phi = v21.Phi() / TMath::Pi() * 180.;
1164  Double_t h1 = bf / 2;
1165  Double_t bl1 = cf / 2;
1166  Double_t tl1 = af / 2;
1167  Double_t alpha1 = -TMath::ATan((cf-af)/2/bf) / TMath::Pi() * 180.;
1168  Double_t h2 = br / 2;
1169  Double_t bl2 = cr / 2;
1170  Double_t tl2 = ar / 2;
1171  Double_t alpha2 = -TMath::ATan((cr-ar)/2/br) / TMath::Pi() * 180.;
1172 
1173  TGeoTrap* trap;
1174  if (type > 0)
1175  trap = new TGeoTrap(name, dz, theta, phi, h1, bl1, tl1, alpha1,
1176  h2, bl2, tl2, alpha2);
1177  else
1178  trap = new TGeoTrap(name, dz, theta, 180-phi, h1, bl1, tl1, -alpha1,
1179  h2, bl2, tl2, -alpha2);
1180  return trap;
1181 }
Int_t i
Definition: run_full.C:25
float Tan(float x)
Definition: PndCAMath.h:165
Double_t
TString name
Double_t crystal_theta[44]
Double_t Pi
void CreateAlveoleShapesAndMatricesZ ( )

Definition at line 775 of file createRootGeoFileBarrel_2018v1.C.

References alveoles_L_minus, alveoles_L_plus, alveoles_R_minus, alveoles_R_plus, BF, BR, CF, CAMath::Cos(), CR, create_trap(), crystal_theta, Double_t, get_trap_vertices(), i, matrices_alveole_minus_left, matrices_alveole_minus_right, matrices_alveole_plus_left, matrices_alveole_plus_right, Pi, safety, CAMath::Sin(), CAMath::Tan(), thickness, X, Y, and Z.

Referenced by geom().

775  {
776  Double_t gap_phi = 0.014;
777  Double_t gap_z = 0.014;
778  Double_t thickness = 0.02;
779  Double_t safety = 1e-4;
780  Double_t dlf = 2.;
781  Double_t dlr = 3.;
782 
783  /* Shapes */
784  TGeoTrap* alveoles_L_out[43];
785  TGeoTrap* alveoles_R_out[43];
786  Int_t i = 0;
787  Double_t top_plane_height[] = {21.773, 20.559, 16.696, 12.582};
788  Double_t top_plane_angle[] = {88.574, 110.433, 131.197, 147.2625};
789  for (Int_t iModule = 0; iModule < 4; iModule++) {
790  for (Int_t iAlveole = 0; iAlveole < 3; iAlveole++) {
791  for (Int_t iCrystal = 0; iCrystal < 4; iCrystal++) {
792  Int_t module_id = iModule + 4;
793  Int_t alveole_id = iModule * 3 + iAlveole + 1;
794  Int_t crystal_id = iCrystal + 1;
795 
796  if (module_id == 1 && iAlveole > 0) break; // SuperModule1 only have 1 alveole
797  if (module_id == 7 && iAlveole > 1) break; // Module7 only have 2 alveole
798  if (module_id == 7 && iAlveole == 1 && iCrystal > 2) break; // The last alveole has 3 crystal colomns
799 
800  Double_t l, dlr_i, af, bf, cf, ar, br, cr;
801  Double_t theta_c = TMath::ATan((CR[0] - CF[0])/20.);
802 
803  // out
804  l = dlf + 20. + dlr;
805  dlr_i = 20. + dlr - (20.*TMath::Cos(crystal_theta[i]) +
806  BF[alveole_id-1]*TMath::Sin(crystal_theta[i]));
807  bf = BF[alveole_id-1]/TMath::Cos(crystal_theta[i]) -
808  (BR[alveole_id-1]*TMath::Sin(TMath::Pi()/2+crystal_theta[i+1]-crystal_theta[i])/TMath::Sin(TMath::Pi()/2-crystal_theta[i+1]) -
809  BF[alveole_id-1]/TMath::Cos(crystal_theta[i])) * dlf /
810  (20.*TMath::Cos(crystal_theta[i])+BF[alveole_id-1]*TMath::Sin(crystal_theta[i])) +
811  (thickness + gap_z - safety) / TMath::Cos(crystal_theta[i+1]) +
812  (thickness + gap_z - safety) / TMath::Cos(crystal_theta[i]);
813  cf = CF[0] - (CR[0] - CF[0]) * dlf / 20.
814  + (thickness + gap_phi - safety) / TMath::Cos(theta_c) + thickness + gap_phi - safety;
815  af = cf;
816  br = BR[alveole_id-1]*TMath::Sin(TMath::Pi()/2+crystal_theta[i+1]-crystal_theta[i])/TMath::Sin(TMath::Pi()/2-crystal_theta[i+1]) +
817  (BR[alveole_id-1]*TMath::Sin(TMath::Pi()/2+crystal_theta[i+1]-crystal_theta[i])/TMath::Sin(TMath::Pi()/2-crystal_theta[i+1]) -
818  BF[alveole_id-1]/TMath::Cos(crystal_theta[i])) * dlr_i /
819  (20.*TMath::Cos(crystal_theta[i])+BF[alveole_id-1]*TMath::Sin(crystal_theta[i])) +
820  (thickness + gap_z - safety) / TMath::Cos(crystal_theta[i+1]) +
821  (thickness + gap_z - safety) / TMath::Cos(crystal_theta[i]);
822  cr = CR[0] + (CR[0] - CF[0]) * dlr / 20.
823  + (thickness + gap_phi - safety) / TMath::Cos(theta_c) + thickness + gap_phi - safety;
824  ar = cr;
825 
826  TGeoTrap* alveole_L_out = create_trap(TString::Format("alveole%d_L_out", i+1),
827  1, af, bf, cf, ar, br, cr, l, i);
828  TGeoTrap* alveole_R_out = create_trap(TString::Format("alveole%d_R_out", i+1),
829  -1, af, bf, cf, ar, br, cr, l, i);
830  alveoles_L_out[i] = alveole_L_out;
831  alveoles_R_out[i] = alveole_R_out;
832 
833  // top cut plane
834  Double_t height, theta_top_plane;
835  if (iAlveole == 0 && iCrystal == 0) {
836  height = top_plane_height[iModule];
837  theta_top_plane = 90.+crystal_theta[i]/TMath::Pi()*180.-top_plane_angle[iModule];
838  }
839  TVector3 vertices[8];
840  get_trap_vertices(alveole_L_out, vertices);
841 
842  Double_t y_c1 = vertices[0].Y() +
843  (vertices[4].Y()-vertices[0].Y())/l*(dlf+height);
844  Double_t z_c = vertices[0].Z() + height + dlf;
845  TVector3 h0(0., y_c1, z_c);
846  TVector3 c2center(0, 0, 100.);
847  c2center.RotateX(-theta_top_plane/180*TMath::Pi());
848  TVector3 center_c = h0 + c2center;
849  TGeoRotation* rot_top_plane = new TGeoRotation(
850  TString::Format("rot_top_plane_%d", i+1), 90, 0, 90+theta_top_plane, 90, theta_top_plane, 90);
851  TGeoCombiTrans* cmb_top_plane = new TGeoCombiTrans(TString::Format("cmb_top_plane_%d", i+1),
852  center_c.X(), center_c.Y(), center_c.Z(), rot_top_plane);
853  cmb_top_plane->RegisterYourself();
854 
855  TGeoBBox* top_box = new TGeoBBox(
856  TString::Format("cut_box_%d", i+1), 100., 100., 100.);
857 
858  height = height - (bf+(br-bf)/l*(height+dlf))/
859  (1/TMath::Tan(theta_top_plane/180*TMath::Pi())+TMath::Tan(crystal_theta[i+1]));
860 
861  // bottom cut plane
862  // left-
863  get_trap_vertices(alveole_L_out, vertices);
864  h0.SetXYZ(vertices[2].X(), vertices[2].Y(), vertices[2].Z() + dlf - 0.158);
865  c2center.SetXYZ(0., 0., -100.);
866  c2center.RotateY(-4./180*TMath::Pi());
867  center_c = h0 + c2center;
868  TGeoRotation* rot_bottom_plane_left_minus = new TGeoRotation(
869  TString::Format("rot_bottom_plane_left_minus_%d", i+1), 90-4, 0, 90, 90, 4, 180);
870  TGeoCombiTrans* cmb_bottom_plane_left_minus = new TGeoCombiTrans(TString::Format("cmb_bottom_plane_left_minus_%d", i+1),
871  center_c.X(), center_c.Y(), center_c.Z(), rot_bottom_plane_left_minus);
872  cmb_bottom_plane_left_minus->RegisterYourself();
873 
874  // right-
875  get_trap_vertices(alveole_R_out, vertices);
876  h0.SetXYZ(vertices[2].X(), vertices[2].Y(), vertices[2].Z() + dlf - 0.158);
877  c2center.SetXYZ(0., 0., -100.);
878  c2center.RotateY(-4./180*TMath::Pi());
879  center_c = h0 + c2center;
880  TGeoRotation* rot_bottom_plane_right_minus = new TGeoRotation(
881  TString::Format("rot_bottom_plane_right_minus_%d", i+1), 90-4, 0, 90, 90, 4, 180);
882  TGeoCombiTrans* cmb_bottom_plane_right_minus = new TGeoCombiTrans(TString::Format("cmb_bottom_plane_right_minus_%d", i+1),
883  center_c.X(), center_c.Y(), center_c.Z(), rot_bottom_plane_right_minus);
884  cmb_bottom_plane_right_minus->RegisterYourself();
885 
886  // left+
887  get_trap_vertices(alveole_L_out, vertices);
888  h0.SetXYZ(vertices[1].X(), vertices[1].Y(), vertices[1].Z() + dlf - 0.158);
889  c2center.SetXYZ(0., 0., -100.);
890  c2center.RotateY(4./180*TMath::Pi());
891  center_c = h0 + c2center;
892  TGeoRotation* rot_bottom_plane_left_plus = new TGeoRotation(
893  TString::Format("rot_bottom_plane_left_plus_%d", i+1), 90+4, 0, 90, 90, 4, 0);
894  TGeoCombiTrans* cmb_bottom_plane_left_plus = new TGeoCombiTrans(TString::Format("cmb_bottom_plane_left_plus_%d", i+1),
895  center_c.X(), center_c.Y(), center_c.Z(), rot_bottom_plane_left_plus);
896  cmb_bottom_plane_left_plus->RegisterYourself();
897 
898  // right+
899  get_trap_vertices(alveole_R_out, vertices);
900  h0.SetXYZ(vertices[1].X(), vertices[1].Y(), vertices[1].Z() + dlf - 0.158);
901  c2center.SetXYZ(0., 0., -100.);
902  c2center.RotateY(4./180*TMath::Pi());
903  center_c = h0 + c2center;
904  TGeoRotation* rot_bottom_plane_right_plus = new TGeoRotation(
905  TString::Format("rot_bottom_plane_right_plus_%d", i+1), 90+4, 0, 90, 90, 4, 0);
906  TGeoCombiTrans* cmb_bottom_plane_right_plus = new TGeoCombiTrans(TString::Format("cmb_bottom_plane_right_plus_%d", i+1),
907  center_c.X(), center_c.Y(), center_c.Z(), rot_bottom_plane_right_plus);
908  cmb_bottom_plane_right_plus->RegisterYourself();
909 
910  // in
911  af -= thickness * 2;
912  bf -= thickness * 2;
913  cf -= thickness * 2;
914  ar -= thickness * 2;
915  br -= thickness * 2;
916  cr -= thickness * 2;
917  l += 0.01;
918 
919  TGeoTrap* alveole_L_in = create_trap(TString::Format("alveole%d_L_in", i+1),
920  1, af, bf, cf, ar, br, cr, l, i);
921  TGeoTrap* alveole_R_in = create_trap(TString::Format("alveole%d_R_in", i+1),
922  -1, af, bf, cf, ar, br, cr, l, i);
923 
924  // composite shapes
925  if (i < 28) {
926  alveoles_L_minus[i] = new TGeoCompositeShape(TString::Format("alveole%d_L", i+1),
927  TString::Format("alveole%d_L_out - alveole%d_L_in - cut_box_%d:cmb_top_plane_%d - cut_box_%d:cmb_bottom_plane_left_minus_%d", i+1, i+1, i+1, i+1, i+1, i+1));
928  alveoles_R_minus[i] = new TGeoCompositeShape(TString::Format("alveole%d_R", i+1),
929  TString::Format("alveole%d_R_out - alveole%d_R_in - cut_box_%d:cmb_top_plane_%d - cut_box_%d:cmb_bottom_plane_right_minus_%d", i+1, i+1, i+1, i+1, i+1, i+1));
930  }
931  alveoles_L_plus[i] = new TGeoCompositeShape(TString::Format("alveole%d_L", i+1),
932  TString::Format("alveole%d_L_out - alveole%d_L_in - cut_box_%d:cmb_top_plane_%d - cut_box_%d:cmb_bottom_plane_left_plus_%d", i+1, i+1, i+1, i+1, i+1, i+1));
933  alveoles_R_plus[i] = new TGeoCompositeShape(TString::Format("alveole%d_R", i+1),
934  TString::Format("alveole%d_R_out - alveole%d_R_in - cut_box_%d:cmb_top_plane_%d - cut_box_%d:cmb_bottom_plane_right_plus_%d", i+1, i+1, i+1, i+1, i+1, i+1));
935 
936  i++;
937  } // crystal
938  } // module (alveole)
939  } // supermodule
940 
941  /* Matrices */
942  for (Int_t type = 0; type < 4; type++) { // type: 0 left-, 1 right+, 2 right-, 3 left+
943  const Double_t crystal_gap = 2*safety;
944  const Double_t alveole_gap = 0.09 - 2*thickness - 2*gap_z + 2*safety;
945  Double_t module_gap[] = {
946  0.12 - thickness - gap_z + safety,
947  0.25 - 2*thickness - 2*gap_z + 2*safety,
948  0.3 - 2*thickness - 2*gap_z + 2*safety,
949  0.33 - 2*thickness - 2*gap_z + 2*safety
950  }; // half of the first gap
951 
952  Double_t trans_z = 0;
953 
954  Int_t nModule;
955  if (type == 0 || type == 2) nModule = 3;
956  else nModule = 4;
957  Int_t i = 0;
958  for (Int_t iModule = 0; iModule < nModule; iModule++) {
959  for (Int_t iAlveole = 0; iAlveole < 3; iAlveole++) {
960  for (Int_t iCrystal = 0; iCrystal < 4; iCrystal++) {
961  Int_t module_id;
962  if (type == 0 || type == 2) module_id = 3 - iModule;
963  else module_id = iModule + 4;
964  Int_t alveole_id = iModule * 3 + iAlveole + 1;
965  Int_t crystal_id = iCrystal + 1;
966 
967  if (module_id == 1 && iAlveole > 0) break; // SuperModule1 only have 1 alveole
968  if (module_id == 7 && iAlveole > 1) break; // Module7 only have 2 alveole
969  if (module_id == 7 && iAlveole == 1 && iCrystal > 2) break; // The last alveole has 3 crystal colomns
970 
971  // calculate crystal type
972  TGeoTrap* crystal;
973  if (type == 0 || type == 3) crystal = alveoles_L_out[i];
974  else crystal = alveoles_R_out[i];
975 
976  Double_t gap = crystal_gap;
977  if (iCrystal == 0) {
978  gap = alveole_gap;
979  if (iAlveole == 0) gap = module_gap[iModule];
980  }
981 
982  // calculate the center of the crystal
983  TVector3 pivot2center;
984  TVector3 vertices[8];
985  get_trap_vertices(crystal, vertices);
986  Double_t bf = crystal->GetH1()*2;
987  if (type == 0) {
988  trans_z -= gap/TMath::Cos(crystal_theta[i]) + bf; // trans_z < 0
989  pivot2center = -vertices[1]; // rotate c2p from ROOT coordinate to experimental coordinate
990  }
991  else if (type == 1) {
992  trans_z += gap/TMath::Cos(crystal_theta[i]) + bf; // trans_z < 0
993  pivot2center = -vertices[2]; // rotate c2p from ROOT coordinate to experimental coordinate
994  }
995  else if (type == 2) {
996  trans_z -= gap/TMath::Cos(crystal_theta[i]) + bf; // trans_z < 0
997  pivot2center = -vertices[2];
998  }
999  else {
1000  trans_z += gap/TMath::Cos(crystal_theta[i]) + bf; // trans_z < 0
1001  pivot2center = -vertices[1];
1002  }
1003  if (type == 0) {
1004  pivot2center.RotateX(-TMath::Pi()/2.);
1005  pivot2center.RotateZ(-TMath::Pi()/2.);
1006  }
1007  else if (type == 1) {
1008  pivot2center.RotateZ(TMath::Pi());
1009  pivot2center.RotateX(-TMath::Pi()/2.);
1010  pivot2center.RotateZ(-TMath::Pi()/2.);
1011  }
1012  else if (type == 2) {
1013  pivot2center.RotateZ(TMath::Pi());
1014  pivot2center.RotateX(TMath::Pi()/2.);
1015  pivot2center.RotateZ(TMath::Pi()/2.);
1016  }
1017  else {
1018  pivot2center.RotateX(TMath::Pi()/2.);
1019  pivot2center.RotateZ(TMath::Pi()/2.);
1020  }
1021  TVector3 pivot;
1022  if (type == 0 || type == 1)
1023  pivot = TVector3(57.+0.112126393-dlf,
1024  -2*safety-(gap_phi+thickness-safety),
1025  trans_z+3.7); // EMC radius is 570mm, crystal offset w.r.t. target is 37mm
1026  else
1027  pivot = TVector3(57.-dlf, -(gap_phi+thickness-safety),
1028  trans_z+3.7);
1029  TVector3 center = pivot + pivot2center;
1030 
1031  // calculate the rotation matrix and place the crystal
1032  TGeoRotation* crystal_rot;
1033  char avleole_sign;
1034  char crystal_type;
1035  if (type == 0) {
1036  avleole_sign = 'm';
1037  crystal_type = 'L';
1038  }
1039  else if (type == 1) {
1040  avleole_sign = 'p';
1041  crystal_type = 'R';
1042  }
1043  else if (type == 1) {
1044  avleole_sign = 'm';
1045  crystal_type = 'R';
1046  }
1047  else {
1048  avleole_sign = 'p';
1049  crystal_type = 'L';
1050  }
1051  if (type == 0) {
1052  crystal_rot = new TGeoRotation(
1053  TString::Format("rot-alveole-%d%c-%c%d",
1054  alveole_id, avleole_sign, crystal_type, crystal_id),
1055  90, 270, 180, 180, 90, 0);
1056  }
1057  else if (type == 1) {
1058  crystal_rot = new TGeoRotation(
1059  TString::Format("rot-alveole-%d%c-%c%d",
1060  alveole_id, avleole_sign, crystal_type, crystal_id),
1061  90, 90, 0, 0, 90, 0);
1062  }
1063  else if (type == 2) {
1064  crystal_rot = new TGeoRotation(
1065  TString::Format("rot-alveole-%d%c-%c%d",
1066  alveole_id, avleole_sign, crystal_type, crystal_id),
1067  90, 270, 180, 180, 90, 0);
1068  }
1069  else {
1070  crystal_rot = new TGeoRotation(
1071  TString::Format("rot-alveole-%d%c-%c%d",
1072  alveole_id, avleole_sign, crystal_type, crystal_id),
1073  90, 90, 0, 0, 90, 0);
1074  }
1075  crystal_rot->RegisterYourself();
1076 
1077 
1078  TGeoCombiTrans* crystal_ctrans =
1079  new TGeoCombiTrans(TString::Format("trans-alveole-%d%c-%c%d",
1080  alveole_id, avleole_sign, crystal_type, crystal_id),
1081  center.X(), center.Y(), center.Z(), crystal_rot);
1082 
1083  if (type == 0)
1084  matrices_alveole_minus_left[i] = crystal_ctrans;
1085  else if (type == 1)
1086  matrices_alveole_plus_right[i] = crystal_ctrans;
1087  else if (type == 2)
1088  matrices_alveole_minus_right[i] = crystal_ctrans;
1089  else
1090  matrices_alveole_plus_left[i] = crystal_ctrans;
1091 
1092  i++;
1093  } // crystal
1094  } // module (alveole)
1095  } // supermodule
1096  } // type
1097 }
TGeoMatrix * matrices_alveole_minus_right[28]
Double_t CR[11]
void get_trap_vertices(TGeoTrap *, TVector3 *points)
TGeoMatrix * matrices_alveole_plus_right[43]
Int_t i
Definition: run_full.C:25
Double_t BF[11]
TGeoMatrix * matrices_alveole_plus_left[43]
static T Sin(const T &x)
Definition: PndCAMath.h:42
TGeoTrap * create_trap(TString name, Int_t type, Double_t af, Double_t bf, Double_t cf, Double_t ar, Double_t br, Double_t cr, Double_t l)
float Tan(float x)
Definition: PndCAMath.h:165
TGeoMatrix * matrices_alveole_minus_left[28]
Double_t thickness
double Y
Definition: anaLmdDigi.C:68
TGeoCompositeShape * alveoles_R_plus[43]
static T Cos(const T &x)
Definition: PndCAMath.h:43
Double_t
TGeoCompositeShape * alveoles_L_minus[28]
Double_t BR[11]
TGeoCompositeShape * alveoles_R_minus[28]
double X
Definition: anaLmdDigi.C:68
#define safety
Definition: createSTT.C:44
double Z
Definition: anaLmdDigi.C:68
TGeoCompositeShape * alveoles_L_plus[43]
Double_t crystal_theta[44]
Double_t Pi
Double_t CF[11]
void CreateFrontInsertShapesAndMatricesZ ( )

Definition at line 698 of file createRootGeoFileBarrel_2018v1.C.

References AF, AR, BF, BR, CF, CR, create_trap(), crystals_L, crystals_R, Double_t, front_inserts_L, front_inserts_R, get_trap_vertices(), i, matrices_fe2cry_left, matrices_fe2cry_right, and thickness.

Referenced by geom().

698  {
699  // shapes
700  TGeoTrap* front_inserts_L_out[11];
701  TGeoTrap* front_inserts_R_out[11];
702  for (Int_t i = 0; i < 11; i++) {
703  Double_t l, dl, af, bf, cf, ar, br, cr, thickness;
704 
705  dl = 0.0065 + 0.00001; // wrapping + safe gap
706  cr = CF[i] - (CR[i] - CF[i]) * dl / 20.;
707  //ar = cr;
708  ar = AF[i] - (AR[i] - AF[i]) * dl / 20.;
709  br = BF[i] - (BR[i] - BF[i]) * dl / 20.;
710 
711  dl += 0.3;
712  cf = CF[i] - (CR[i] - CF[i]) * dl / 20.;
713  //af = cf;
714  af = AF[i] - (AR[i] - AF[i]) * dl / 20.;
715  bf = BF[i] - (BR[i] - BF[i]) * dl / 20.;
716 
717  l = 0.3;
718 
719  front_inserts_L_out[i] = create_trap(TString::Format("front_insert_%d_L_out", i+1),
720  1, af, bf, cf, ar, br, cr, l);
721  front_inserts_R_out[i] = create_trap(TString::Format("front_insert_%d_R_out", i+1),
722  -1, af, bf, cf, ar, br, cr, l);
723 
724  TGeoTube* big_hole = new TGeoTube("front_insert_big_hole", 0., 0.6, 0.151);
725  TGeoTube* small_hole = new TGeoTube("front_insert_small_hole", 0., 0.2, 0.151);
726  TGeoTranslation* trans_big_hole_L = new TGeoTranslation("trans_front_insert_big_hole_L", -0.13, 0.13, 0.);
727  trans_big_hole_L->RegisterYourself();
728  TGeoTranslation* trans_small_hole_L = new TGeoTranslation("trans_front_insert_small_hole_L", 0.67, -0.67, 0.);
729  trans_small_hole_L->RegisterYourself();
730  TGeoTranslation* trans_big_hole_R = new TGeoTranslation("trans_front_insert_big_hole_R", 0.13, 0.13, 0.);
731  trans_big_hole_R->RegisterYourself();
732  TGeoTranslation* trans_small_hole_R = new TGeoTranslation("trans_front_insert_small_hole_R", -0.67, -0.67, 0.);
733  trans_small_hole_R->RegisterYourself();
734 
735  front_inserts_L[i] = new TGeoCompositeShape(TString::Format("front_insert_%d_L", i+1),
736  TString::Format("front_insert_%d_L_out - front_insert_big_hole:trans_front_insert_big_hole_L", i+1));
737  //TString::Format("front_insert_%d_L_out - front_insert_big_hole:trans_front_insert_big_hole_L - front_insert_small_hole:trans_front_insert_small_hole_L", i+1));
738  front_inserts_R[i] = new TGeoCompositeShape(TString::Format("front_insert_%d_R", i+1),
739  TString::Format("front_insert_%d_R_out - front_insert_big_hole:trans_front_insert_big_hole_R", i+1));
740  //TString::Format("front_insert_%d_R_out - front_insert_big_hole:trans_front_insert_big_hole_R - front_insert_small_hole:trans_front_insert_small_hole_R", i+1));
741 
742  }
743 
744  // matrices
745  for (Int_t i = 0; i < 11; i++) {
746  // left
747  TGeoTrap* crystal = crystals_L[i];
748  TGeoTrap* front_insert = front_inserts_L_out[i];
749 
750  TVector3 vertices[8];
751  get_trap_vertices(crystal, vertices);
752  TVector3 ccry_to_v0cry = vertices[0]; // crystal's center point to vertex 0
753  TVector3 v0cry_to_v4ins = TVector3(0., 0., -(0.0065 + 0.00001)); // crystal's vertex 0 to insert's vertex 4
754  get_trap_vertices(front_insert, vertices);
755  TVector3 v4ins_to_cins = -vertices[4]; // insert's vertex 4 to center point
756  TVector3 ccry_to_cins = ccry_to_v0cry + v0cry_to_v4ins + v4ins_to_cins; // crystal's certer to insert's center
757 
758  matrices_fe2cry_left[i] = new TGeoTranslation(ccry_to_cins.X(), ccry_to_cins.Y(), ccry_to_cins.Z());
759 
760  // right
761  crystal = crystals_R[i];
762  front_insert = front_inserts_R_out[i];
763 
764  get_trap_vertices(crystal, vertices);
765  TVector3 ccry_to_v3cry = vertices[3]; // crystal's center point to vertex 3
766  TVector3 v3cry_to_v7ins = TVector3(0., 0., -(0.0065 + 0.00001)); // crystal's vertex 3 to insert's vertex 7
767  get_trap_vertices(front_insert, vertices);
768  TVector3 v7ins_to_cins = -vertices[7]; // insert's vertex 7 to center point
769  ccry_to_cins = ccry_to_v3cry + v3cry_to_v7ins + v7ins_to_cins; // crystal's certer to insert's center
770 
771  matrices_fe2cry_right[i] = new TGeoTranslation(ccry_to_cins.X(), ccry_to_cins.Y(), ccry_to_cins.Z());
772  }
773 }
TGeoTranslation * matrices_fe2cry_left[11]
Double_t CR[11]
void get_trap_vertices(TGeoTrap *, TVector3 *points)
TGeoTrap * crystals_L[11]
Int_t i
Definition: run_full.C:25
Double_t BF[11]
TGeoTrap * create_trap(TString name, Int_t type, Double_t af, Double_t bf, Double_t cf, Double_t ar, Double_t br, Double_t cr, Double_t l)
Double_t thickness
Double_t AR[11]
Double_t AF[11]
Double_t
TGeoCompositeShape * front_inserts_L[11]
Double_t BR[11]
TGeoCompositeShape * front_inserts_R[11]
TGeoTranslation * matrices_fe2cry_right[11]
TGeoTrap * crystals_R[11]
Double_t CF[11]
void geom ( )

Definition at line 79 of file createRootGeoFileBarrel_2018v1.C.

References ABS, AF, AR, BF, BR, CalculateCrystalMatricesPhi(), CalculateCrystalMatricesZ(), CF, ConstructSlice(), ConstructTargetSlice(), CR, create_trap(), CreateAlveoleShapesAndMatricesZ(), CreateFrontInsertShapesAndMatricesZ(), crystals_L, crystals_R, Double_t, f, geobuild, geoFace, geoLoad, geoMan, i, Media, rot, theta, thickness, top, wrappings_L, and wrappings_R.

Referenced by FairEvtFilterOnSingleParticleCounts::AcceptGeometry(), FairEvtFilterOnSingleParticleCounts::AndMinMaxGeom(), CrystalShapeTest(), hypGe_DoubleGeo(), hypGe_GeoBuilder_template(), hypGe_TripleGeo(), hypGeCableAbsorptionTest(), hypGeGeantTestGeometry(), hypGeGeoBuilderDEGASBall40_6sym_offset20(), hypGeGeoBuilderDEGASStraight40_offset15(), hypGeGeoBuilderDEGASStraight40_offset20(), hypGeGeoBuilderDouble30cmRadius(), hypGeGeoBuilderDouble30cmRadius_test(), hypGeGeoBuilderSingle(), hypGeGeoBuilderTriple30cmRadius(), hypGeGeoBuilderTriple30cmRadius_test(), hypGeGeoBuilderTriple30cmRadiusCrystalsOnly(), hypGeGeoBuilderTripleBall40Offset10Geometry(), hypGeGeoBuilderTripleBall40Offset10Geometry_STTFitting(), hypGeGeoBuilderTripleBall40Offset10GeometryCrystalsOnly(), hypGeGeoBuilderTripleBall40Offset20Geometry(), hypGeGeoBuilderTripleBall40Offset20Geometry_STTFitting(), hypGeGeoBuilderTripleBall40Offset20Geometry_STTFittingCrystalsOnly(), hypGeGeoBuilderTripleBall40Offset20GeometryCrystalsOnly(), hypGeGeoBuilderTripleStraightGeometry(), hypGeGeoBuilderTripleStraightGeometryCrystalsOnly(), hypGeGeoCOSYBeamDumpTOF(), hypGeGeoCOSYGermaniums(), hypGeGeoCOSYInBeamStuff(), hypGeGeoCOSYsetup2014(), hypGeGeoCOSYsetup2014Actives(), hypGeGeoCOSYsetup2014Passives(), hypGeGeoCOSYSiPm(), hypGeGeoCOSYTarget(), main(), muon_barrel_strip_5bis(), sim(), and simLut().

79  {
80  //--- Definition of a simple geometry
81  gSystem->Load("libGeom");
82  gSystem->Load("libGeoBase");
83  gSystem->Load("libParBase");
84  gSystem->Load("libBase");
85  gSystem->Load("libPndData");
86  gSystem->Load("libPassive");
87 
88  // load create materials/media
89  FairGeoLoader* geoLoad = new FairGeoLoader("TGeo", "FairGeoLoader");
90  FairGeoInterface* geoFace = geoLoad->getGeoInterface();
91  geoFace->setMediaFile("media_pnd.geo");
92  geoFace->readMedia();
93  geoFace->print();
94 
95  FairGeoMedia* Media = geoFace->getMedia();
96  FairGeoBuilder* geobuild = geoLoad->getGeoBuilder();
97 
98  FairGeoMedium* PWO = Media->getMedium("PWO");
99  FairGeoMedium* Prepreg = Media->getMedium("Prepreg");
100  FairGeoMedium* VM2000 = Media->getMedium("VM2000");
101  FairGeoMedium* ABS = Media->getMedium("ABS");
102  FairGeoMedium* Polyurethane = Media->getMedium("Polyurethane");
103  FairGeoMedium* Aluminum = Media->getMedium("Aluminum");
104  geobuild->createMedium(PWO);
105  geobuild->createMedium(Prepreg);
106  geobuild->createMedium(VM2000);
107  geobuild->createMedium(ABS);
108  geobuild->createMedium(Polyurethane);
109  geobuild->createMedium(Aluminum);
110 
111  //--- make the top container volume
112  geoMan = (TGeoManager*)gROOT->FindObject("FAIRGeom");
113  TGeoVolume *top = new TGeoVolumeAssembly("BarrelEMC");
114  geoMan->SetTopVolume(top);
115 
116  /* Crystals */
117  for (Int_t i = 0; i < 11; i++) {
118  crystals_L[i] = create_trap(TString::Format("crystal%d_L", i+1), 1,
119  AF[i], BF[i], CF[i], AR[i], BR[i], CR[i], 20.);
120  crystals_R[i] = create_trap(TString::Format("crystal%d_R", i+1), -1,
121  AF[i], BF[i], CF[i], AR[i], BR[i], CR[i], 20.);
122  }
125 
126  /* Wrapping */
127  for (Int_t i = 0; i < 11; i++) {
128  Double_t l, dl, af, bf, cf, ar, br, cr, thickness;
129 
130  // out
131  thickness = 0.0065;
132  l = 20. + thickness*2;
133  dl = (l - 20.) / 2;
134  af = AF[i] - (AR[i] - AF[i]) * dl / 20. + thickness * 2;
135  bf = BF[i] - (BR[i] - BF[i]) * dl / 20. + thickness * 2;
136  cf = CF[i] - (CR[i] - CF[i]) * dl / 20. + thickness * 2;
137  ar = AR[i] + (AR[i] - AF[i]) * dl / 20. + thickness * 2;
138  br = BR[i] + (BR[i] - BF[i]) * dl / 20. + thickness * 2;
139  cr = CR[i] + (CR[i] - CF[i]) * dl / 20. + thickness * 2;
140 
141  TGeoTrap* wrapping_L_out = create_trap(TString::Format("wrapping%d_L_out", i+1),
142  1, af, bf, cf, ar, br, cr, l);
143  TGeoTrap* wrapping_R_out = create_trap(TString::Format("wrapping%d_R_out", i+1),
144  -1, af, bf, cf, ar, br, cr, l);
145 
146  // in
147  thickness = 0.00001;
148  l = 20. + thickness*2;
149  dl = (l - 20.) / 2;
150  af = AF[i] - (AR[i] - AF[i]) * dl / 20. + thickness * 2;
151  bf = BF[i] - (BR[i] - BF[i]) * dl / 20. + thickness * 2;
152  cf = CF[i] - (CR[i] - CF[i]) * dl / 20. + thickness * 2;
153  ar = AR[i] + (AR[i] - AF[i]) * dl / 20. + thickness * 2;
154  br = BR[i] + (BR[i] - BF[i]) * dl / 20. + thickness * 2;
155  cr = CR[i] + (CR[i] - CF[i]) * dl / 20. + thickness * 2;
156 
157  TGeoTrap* wrapping_L_in = create_trap(TString::Format("wrapping%d_L_in", i+1),
158  1, af, bf, cf, ar, br, cr, l);
159  TGeoTrap* wrapping_R_in = create_trap(TString::Format("wrapping%d_R_in", i+1),
160  -1, af, bf, cf, ar, br, cr, l);
161 
162  // boolean shape
163  wrappings_L[i] = new TGeoCompositeShape(TString::Format("wrapping%d_L", i+1),
164  TString::Format("wrapping%d_L_out - wrapping%d_L_in", i+1, i+1));
165  wrappings_R[i] = new TGeoCompositeShape(TString::Format("wrapping%d_R", i+1),
166  TString::Format("wrapping%d_R_out - wrapping%d_R_in", i+1, i+1));
167  }
168 
169  /* Front inserts */
171 
172  /* Alveoles */
174 
175  /* construct whole barrel */
176  TGeoVolumeAssembly* Slice = new TGeoVolumeAssembly("Slice");
177  ConstructSlice(Slice);
178  TGeoVolumeAssembly* Slice_target = new TGeoVolumeAssembly("Slice_target");
179  ConstructTargetSlice(Slice_target);
180  for (Int_t i = 0; i < 1; i++) {
181  double theta = 360./16 * i;
182  TGeoRotation* rot = new TGeoRotation(TString::Format("rot_slice_%d", i),
183  90, theta, 90, 90+theta, 0, 0);
184  if (i == 0 || i == 8) top->AddNode(Slice_target, i+1, rot);
185  else top->AddNode(Slice, i+1, rot);
186  }
187 
188  //--- close the geometry
189  geoMan->CloseGeometry();
190 
191  /* draw geometry */
192  geoMan->SetVisLevel(4);
193  geoMan->SetVisOption(0);
194  top->Draw("ogl");
195  //TView* view = gPad->GetView();
196  //view->ShowAxis();
197 
198  /* overlap checking */
199  //geoMan->CheckOverlaps(1e-5);
200 
201  /* raytrace */
202  //raytrace_z(geoMan, "emc_raytrace_xz.root");
203  //raytrace_x(geoMan, "emc_raytrace_xy.root");
204 
205  TFile* f = new TFile("emc_module12_2018v1.root", "recreate");
206  top->Write();
207  f->Close();
208 }
Double_t CR[11]
TGeoTrap * crystals_L[11]
FairGeoLoader * geoLoad
Int_t i
Definition: run_full.C:25
FairGeoMedia * Media
Double_t BF[11]
void CalculateCrystalMatricesZ()
TGeoManager * geoMan
TGeoCompositeShape * wrappings_R[11]
#define ABS
TGeoTrap * create_trap(TString name, Int_t type, Double_t af, Double_t bf, Double_t cf, Double_t ar, Double_t br, Double_t cr, Double_t l)
Double_t thickness
void CreateAlveoleShapesAndMatricesZ()
Double_t AR[11]
TGeoVolume * top
void ConstructSlice(TGeoVolume *vol)
TGeoCompositeShape * wrappings_L[11]
void CreateFrontInsertShapesAndMatricesZ()
FairGeoBuilder * geobuild
Double_t AF[11]
Double_t
Double_t BR[11]
TFile * f
Definition: bump_analys.C:12
void CalculateCrystalMatricesPhi()
void ConstructTargetSlice(TGeoVolume *vol)
TGeoRotation rot
FairGeoInterface * geoFace
TGeoTrap * crystals_R[11]
Double_t CF[11]
void get_trap_vertices ( TGeoTrap *  trap,
TVector3 *  points 
)

Definition at line 1183 of file createRootGeoFileBarrel_2018v1.C.

References Double_t, i, and p.

Referenced by CalculateCrystalMatricesZ(), CreateAlveoleShapesAndMatricesZ(), and CreateFrontInsertShapesAndMatricesZ().

1183  {
1184  Double_t xy[8][2];
1185  Double_t* p = trap->GetVertices();
1186  Int_t off = 0;
1187  for (Int_t i = 0; i < 8; i++) {
1188  for (Int_t j = 0; j < 2; j++) {
1189  xy[i][j] = *(p + off);
1190  off++;
1191  }
1192  }
1193  for (Int_t i = 0; i < 8; i++) {
1194  if (i < 4) {
1195  points[i].SetXYZ(xy[i][0], xy[i][1], -trap->GetDz());
1196  }
1197  else {
1198  points[i].SetXYZ(xy[i][0], xy[i][1], trap->GetDz());
1199  }
1200  }
1201 }
Double_t p
Definition: anasim.C:58
Int_t i
Definition: run_full.C:25
Double_t
int main ( void  )

Definition at line 75 of file createRootGeoFileBarrel_2018v1.C.

References geom().

75  {
76  geom();
77 }
void raytrace_x ( TGeoManager *  geom,
TString  filename 
)

Definition at line 1247 of file createRootGeoFileBarrel_2018v1.C.

References count, Double_t, dy, dz, f, n, phi, Pi, t, theta, xmin, y, and z.

1247  {
1248  Double_t zmin = 4;
1249  Double_t zmax = 4.5;
1250  Double_t nz = 1E2;
1251  Double_t dz = (zmax - zmin) / nz;
1252  Double_t ymin = -100.;
1253  Double_t ymax = 100.;
1254  Double_t ny = 1E4;
1255  Double_t dy = (ymax - ymin) / ny;
1256  Double_t xmin = 0.;
1257  Double_t theta = 0.1;
1258  Double_t phi = 0.;
1259 
1260  TVector3 n(1., 0., 0.);
1261  n.RotateZ(theta/180.*TMath::Pi());
1262  n.RotateX(phi/180.*TMath::Pi());
1263  TFile* f = new TFile(filename, "recreate");
1264  TTree* t = new TTree("raytrace", "raytrace");
1265  Double_t _x, _y, _z;
1266  t->Branch("x", &_x, "x/D");
1267  t->Branch("y", &_y, "y/D");
1268  t->Branch("z", &_z, "z/D");
1269  Int_t count = 0.;
1270  for (Double_t y = ymin; y < ymax; y += dy) {
1271  for (Double_t z = zmin; z < zmax; z += dz) {
1272  Int_t unit = (Int_t)(ny*nz/100);
1273  if (count % unit == 0) {
1274  cout << count / unit << "% proceeded" << endl;
1275  }
1276  geoMan->InitTrack(xmin, y, z, n.X(), n.Y(), n.Z());
1277  do {
1278  geoMan->FindNextBoundaryAndStep();
1279  _x = geoMan->GetCurrentPoint()[0];
1280  _y = geoMan->GetCurrentPoint()[1];
1281  _z = geoMan->GetCurrentPoint()[2];
1282  t->Fill();
1283  } while (_x < 1000 && _x > -1000);
1284  //} while (!geoMan->IsOutside());
1285  count++;
1286  }
1287  }
1288  t->Write();
1289  f->Close();
1290 }
double dy
TGeoManager * geoMan
int n
Double_t
TFile * f
Definition: bump_analys.C:12
Double_t z
TTree * t
Definition: bump_analys.C:13
int count
Double_t xmin
Double_t y
Double_t Pi
const string filename
void raytrace_z ( TGeoManager *  geom,
TString  filename 
)

Definition at line 1203 of file createRootGeoFileBarrel_2018v1.C.

References count, Double_t, dx, dy, f, n, phi, Pi, t, theta, x, xmax, xmin, and y.

1203  {
1204  Double_t xmin = 50.;
1205  Double_t xmax = 100.;
1206  Double_t nx = 1e2;
1207  Double_t dx = (xmax - xmin) / nx;
1208  Double_t ymin = -5.;
1209  Double_t ymax = 5.;
1210  Double_t ny = 100;
1211  Double_t dy = (ymax - ymin) / ny;
1212  Double_t zmin = -100.;
1213  Double_t theta = 1;
1214  Double_t phi = 90.;
1215 
1216  TVector3 n(0., 0., 1.);
1217  n.RotateY(theta/180.*TMath::Pi());
1218  n.RotateZ(phi/180.*TMath::Pi());
1219  TFile* f = new TFile(filename, "recreate");
1220  TTree* t = new TTree("raytrace", "raytrace");
1221  Double_t _x, _y, _z;
1222  t->Branch("x", &_x, "x/D");
1223  t->Branch("y", &_y, "y/D");
1224  t->Branch("z", &_z, "z/D");
1225  Int_t count = 0.;
1226  for (Double_t x = xmin; x < xmax; x += dx) {
1227  for (Double_t y = ymin; y < ymax; y += dy) {
1228  Int_t unit = (Int_t)(nx*ny/100);
1229  if (count % unit == 0) {
1230  cout << count / unit << "% procecced" << endl;
1231  }
1232  geoMan->InitTrack(x, y, zmin, n.X(), n.Y(), n.Z());
1233  while (!geoMan->IsOutside()) {
1234  geoMan->FindNextBoundaryAndStep(10000., kTRUE);
1235  _x = geoMan->GetCurrentPoint()[0];
1236  _y = geoMan->GetCurrentPoint()[1];
1237  _z = geoMan->GetCurrentPoint()[2];
1238  t->Fill();
1239  }
1240  count++;
1241  }
1242  }
1243  t->Write();
1244  f->Close();
1245 }
double dy
TGeoManager * geoMan
Double_t xmax
int n
Double_t
TFile * f
Definition: bump_analys.C:12
double dx
Double_t x
TTree * t
Definition: bump_analys.C:13
int count
Double_t xmin
Double_t y
Double_t Pi
const string filename

Variable Documentation

Double_t AF[11] = {2.121, 2.118, 2.117, 2.117, 2.117, 2.119, 2.122, 2.123, 2.123, 2.125, 2.125}

Definition at line 48 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CreateFrontInsertShapesAndMatricesZ(), and geom().

TGeoCompositeShape* alveoles_L_minus[28]
TGeoCompositeShape* alveoles_L_plus[43]
TGeoCompositeShape* alveoles_R_minus[28]
TGeoCompositeShape* alveoles_R_plus[43]
Double_t AR[11] = {2.904, 2.878, 2.836, 2.790, 2.735, 2.672, 2.623, 2.570, 2.514, 2.470, 2.435}

Definition at line 51 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CreateFrontInsertShapesAndMatricesZ(), and geom().

Double_t BF[11] = {2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128, 2.128}
Double_t BR[11] = {2.875, 2.875, 2.875, 2.722, 2.722, 2.722, 2.547, 2.547, 2.547, 2.442, 2.442}
Double_t CF[11] = {2.127, 2.139, 2.151, 2.160, 2.169, 2.178, 2.186, 2.191, 2.195, 2.200, 2.202}
Double_t CR[11] = {2.912, 2.907, 2.881, 2.845, 2.801, 2.747, 2.699, 2.651, 2.600, 2.556, 2.523}
Double_t crystal_dz[43]

Definition at line 42 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CalculateCrystalMatricesZ().

Double_t crystal_theta[44]
TGeoTrap* crystals_L[11]
TGeoTrap* crystals_R[11]
TGeoCompositeShape* front_inserts_L[11]
TGeoCompositeShape* front_inserts_R[11]
TGeoManager* geoMan
TGeoMatrix* matrices_alveole_minus_left[28]
TGeoMatrix* matrices_alveole_minus_right[28]
TGeoMatrix* matrices_alveole_plus_left[43]
TGeoMatrix* matrices_alveole_plus_right[43]
TGeoTranslation* matrices_crystal2alveole_L[11]

Definition at line 45 of file createRootGeoFileBarrel_2018v1.C.

TGeoTranslation* matrices_crystal2alveole_R[11]

Definition at line 46 of file createRootGeoFileBarrel_2018v1.C.

TGeoTranslation* matrices_fe2cry_left[11]
TGeoTranslation* matrices_fe2cry_right[11]
TGeoCombiTrans* matrices_minus_left[28]

Definition at line 32 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CalculateCrystalMatricesZ(), and ConstructModule().

TGeoCombiTrans* matrices_minus_right[28]

Definition at line 33 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CalculateCrystalMatricesZ(), and ConstructModule().

TGeoMatrix* matrices_phi[5]
TGeoCombiTrans* matrices_plus_left[43]

Definition at line 34 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CalculateCrystalMatricesZ(), and ConstructModule().

TGeoCombiTrans* matrices_plus_right[43]

Definition at line 35 of file createRootGeoFileBarrel_2018v1.C.

Referenced by CalculateCrystalMatricesZ(), and ConstructModule().

Double_t theta_hole_corr = -(4.997 + 3.862213309 + 90.)
TGeoCompositeShape* wrappings_L[11]

Definition at line 24 of file createRootGeoFileBarrel_2018v1.C.

Referenced by ConstructModule(), and geom().

TGeoCompositeShape* wrappings_R[11]

Definition at line 25 of file createRootGeoFileBarrel_2018v1.C.

Referenced by ConstructModule(), and geom().