]> jspc29.x-matter.uni-frankfurt.de Git - mvd_geometry.git/commitdiff
v0.5
authorPhilipp Klaus <klaus@physik.uni-frankfurt.de>
Fri, 4 Aug 2017 15:49:31 +0000 (17:49 +0200)
committerPhilipp Klaus <klaus@physik.uni-frankfurt.de>
Fri, 4 Aug 2017 15:49:31 +0000 (17:49 +0200)
mvd.C

diff --git a/mvd.C b/mvd.C
index 145851fdca19f13aad1b4442ec7c6d0c19862067..c3e6ac320d49cf4c2665ce861716ef50c022f751 100644 (file)
--- a/mvd.C
+++ b/mvd.C
@@ -1,82 +1,96 @@
-int numStations = 4;
+static const int numStations = 4;
 Float_t quadrantBeamOffset[4] = {-0.54, -0.54, -0.82, -1.04};
-Float_t carrierDimensions[4][2] = {{5.0, 5.15}, {8.0, 8.0}, {11.0, 9.9}, {14.0, 12.75}};
+Float_t carrierDimensions[4][3] = {{5.0, 5.15, 0.015}, {8.0, 8.0, 0.015}, {11.0, 9.9, 0.0375}, {14.0, 12.75, 0.0375}};
 Float_t stationPosition[4] = {5.0, 10.0, 15.0, 20.0};
 int sensorRows[4] = {3, 5, 7, 10};
 int sensorCols[4] = {1, 2, 3,  4};
 // 3xFSBB : 30 x 13 mm2 (of which the lower 3mm inactive)
 Float_t sensorDimensionsActive[3]  = {3., 1.0, 0.005};
 Float_t sensorDimensionsPassive[3] = {3., 0.3, 0.005};
-Float_t sensorPitch = sensorDimensionsActive[1] + sensorDimensionsPassive[1] - 0.05;
+Float_t sensorPitch = sensorDimensionsActive[1] - 0.05;
+Float_t sensorSpacing = 0.01;
 Float_t fpcDimensions[3] = {-1, 1.0, 0.004};
+Float_t explosion = .0; // set to 0. for no explosion, and 5. for a heavy explosion
+
+TGeoVolume *top;
+TGeoVolume *stations[numStations];
 
 void mvd()
    {   
    TGeoManager *manager = new TGeoManager("Chamber Layout", "Chamber Layout");
    TGeoMaterial *mat = new TGeoMaterial("Vacuum", 0.,0.,0.);
    TGeoMaterial *tmat = new TGeoMaterial("TubeMat", 26.98,13.,2.700);
-   tmat->SetTransparency(50);
+   tmat->SetTransparency(0);
    TGeoMaterial *wmat = new TGeoMaterial("WireMat", 63.546,29.,8.920);
    TGeoMaterial *matAl = new TGeoMaterial("Al", 26.98,13,2.7);
+   matAl->SetTransparency(0);
    TGeoMedium *med = new TGeoMedium("Vacuum", 0, mat);
    TGeoMedium *tmed = new TGeoMedium("TubeMed", 0, tmat);
    TGeoMedium *wmed = new TGeoMedium("WireMed", 0, wmat);
    TGeoMedium *Al = new TGeoMedium("AlMed", 2, matAl);
-   TGeoVolume *top = manager->MakeBox("Top",med,150.,150.,150.);
+   top = manager->MakeBox("Top",med,300.,300.,300.);
 
    manager->SetTopVolume(top);
 
-   TGeoVolume *sensor_active  = manager->MakeBox("sensor_active",  Al, sensorDimensionsActive[0]/2, sensorDimensionsActive[1]/2, sensorDimensionsActive[2]/2);
-   TGeoVolume *sensor_passive = manager->MakeBox("sensor_passive", Al, sensorDimensionsActive[0]/2, sensorDimensionsActive[1]/2, sensorDimensionsActive[2]/2);
-   sensor_active->SetTransparency(10);
-   sensor_active->SetLineColor(20);
-   sensor_passive->SetTransparency(2);
-   sensor_passive->SetLineColor(24);
+   TGeoVolume *sensorActive  = manager->MakeBox("sensorActive",  Al, sensorDimensionsActive[0]/2, sensorDimensionsActive[1]/2, sensorDimensionsActive[2]/2);
+   TGeoVolume *sensorPassive = manager->MakeBox("sensorPassive", Al, sensorDimensionsPassive[0]/2, sensorDimensionsPassive[1]/2, sensorDimensionsPassive[2]/2);
+   sensorActive->SetTransparency(0);
+   sensorActive->SetLineColor(20);
+   sensorPassive->SetTransparency(0);
+   sensorPassive->SetLineColor(24);
 
    TGeoTranslation T;
    TGeoRotation R;
    TGeoCombiTrans *M;
-   char station_name[30], quadrant_name[30], carrier_name[30], quadrant_name[30];
+   char station_name[30], quadrant_name[30], carrier_name[30];
 
    for (int i=0; i<numStations; i++)
       {
       cout << "Station " << i << endl;
-      sprintf(station_name, "S[%i]",i);
+      sprintf(station_name, "station_S%i", i);
       TGeoVolume *station = new TGeoVolumeAssembly(station_name);
       
-      sprintf(quadrant_name, "S[%i]_Q",i);
+      sprintf(quadrant_name, "quadrant_S%i", i);
       TGeoVolume *quadrant = new TGeoVolumeAssembly(quadrant_name);
 
-      Float_t cd[3] = carrierDimensions[i];
-      TGeoVolume* carrier = manager->MakeBox("carrier", Al, cd[0]/2, cd[1]/2, cd[2]/2);
-      //carrier->GetShape()->SetTransform(new TGeoTranslation(8.0/2,8.0/2,0.05/2));
+      sprintf(carrier_name, "carrier_S%i", i);
+      Float_t cd[3] = {carrierDimensions[i][0], carrierDimensions[i][1], carrierDimensions[i][2]};
+      TGeoVolume* carrier = manager->MakeBox(carrier_name, Al, cd[0]/2, cd[1]/2, cd[2]/2);
       carrier->SetTransparency(1);
       carrier->SetLineColor(27);
-      quadrant->AddNode(carrier, 1);
+      T.SetTranslation(-cd[0]/2, -cd[1]/2, 0);
+      R.SetAngles(0.,0.,0.);
+      M = new TGeoCombiTrans(T,R);
+      quadrant->AddNode(carrier, 1, M);
+      cout << "  # of sensors per quadrant: " << sensorCols[i]*sensorRows[i] << endl;
+      // Sensors
       for (int k=0; k<sensorCols[i]; k++)
          {
          for (int l=0; l<sensorRows[i]; l++)
             {
             // Active
-            cout << "  Sensor " << k*sensorCols[i] + l << endl;
-            Float_t z_offset = (-1.0)**(l+1) * (sensorDimensionsActive[2] + cd[2])/2.0;
-            T.SetTranslation(sensorDimensionsActive[0]*(Float_t)k, sensorPitch*l, z_offset);
+            Float_t x_offset = -sensorDimensionsActive[0]*(0.5 + k);
+            Float_t y_offset = -sensorDimensionsActive[1]/2 - sensorPitch*l;
+            Float_t z_offset = TMath::Power(-1, l+1) * (sensorDimensionsActive[2] + cd[2])*(1+explosion);
+            T.SetTranslation(x_offset, y_offset, z_offset);
             R.SetAngles(0.,0.,0.);
             M = new TGeoCombiTrans(T,R);
-            quadrant->AddNode(sensor_active, 1, M);
+            quadrant->AddNode(sensorActive, 1, M);
             // Passive
-            T.SetTranslation(sensorDimensionsActive[0]*(Float_t)k, sensorDimensionsActive[1]+sensorPitch*l, z_offset);
+            y_offset = -sensorDimensionsPassive[1]/2 - sensorDimensionsActive[1] - sensorPitch*l;
+            T.SetTranslation(x_offset, y_offset, z_offset);
             R.SetAngles(0.,0.,0.);
             M = new TGeoCombiTrans(T,R);
-            quadrant->AddNode(sensor_passive, 1, M);
+            quadrant->AddNode(sensorPassive, 1, M);
             }
          }
       for (int j=0; j<4; j++)
          {
          cout << " Quadrant " << j << endl;
-         T.SetTranslation(-quadrantBeamOffset[i]-cd[0]/2, -quadrantBeamOffset[i]-cd[1]/2, 0.0);
-         //T.SetTranslation(-quadrantBeamOffset[i], -quadrantBeamOffset[i], 0.0);
-         R.SetAngles(0.,0.,90.*j);
+         Float_t x_sign = (j == 0 || j == 3) ? -1. :  1.;
+         Float_t y_sign = (j == 0 || j == 1) ?  1. : -1.;
+         T.SetTranslation(x_sign*quadrantBeamOffset[i], y_sign*quadrantBeamOffset[i], (explosion/10.)*(j-1.5));
+         R.SetAngles(0.,0.,-90.*j);
          M = new TGeoCombiTrans(T,R);
          station->AddNode(quadrant, j, M);
          }
@@ -84,13 +98,18 @@ void mvd()
       R.SetAngles(0.,0.,0.);
       M = new TGeoCombiTrans(T,R);
       top->AddNode(station, 1, M);
+      stations[i] = station;
       }
    manager->CloseGeometry();
-   //TCanvas *c3D = new TCanvas("c3D","MVD Layout",1200,800);
-   //top->Draw("ogl");
+   TCanvas *c3D = new TCanvas("c3D","MVD Layout",1600,1200);
+   top->Draw("ogl");
    //station->Draw("ogl");
    //station->Draw("");
-   top->Draw("");
-   TView *view = gPad->GetView()->ShowAxis();
+   //quadrant->Draw("ogl");
+   //quadrant->Draw("");
+   //quadrant->Raytrace();
+   //top->Draw("");
+   gPad->GetView()->ShowAxis();
+   //top->Raytrace();
    top->Export("mvd.root");
    }