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

diff --git a/mvd.C b/mvd.C
index 5a0812dc9969ddc8eeac197d14aae9e2e1d73898..3b4cccbc2d37d7a223befbbb7bed4c59ea32bb95 100644 (file)
--- a/mvd.C
+++ b/mvd.C
@@ -1,13 +1,14 @@
 static const int numStations = 4;
 // all units are in cm
 Float_t quadrantBeamOffset[4] = {0.54, 0.54, 0.82, 1.04};
-Float_t carrierDimensions[4][3] = {{5.0, 5.15, 0.015}, {8.0, 8.0, 0.015}, {11.0, 9.9, 0.0300}, {14.0, 12.75, 0.0300}};
-//Float_t stationPosition[4] = {35.0, 65.0, 110.0, 170.0};
+Float_t carrierDimensions[4][3] = {{4.0, 4.15, 0.015}, {7.0, 7.0, 0.015}, {10.0, 8.9, 0.0300}, {13.0, 11.75, 0.0300}};
+//Float_t stationPosition[4] = {8.0, 12.0, 16.0, 20.0};
 Float_t stationPosition[4] = {5.0, 10.0, 15.0, 20.0};
-Float_t carrierClampOverlap = 1.0;
+Float_t carrierClampOverlap = 0.0;
 Float_t heatsinkWidth[4] = {23.1, 31.1, 34.34, 37.7};
 Float_t heatsinkHeight = 37.7;
 Float_t heatsinkThickness = 1.;
+Float_t fpcWidth = 1.9;
 Float_t fpcThickness = 0.0009;
 Float_t glueThickness = 0.0008;
 int sensorRows[4] = {2, 5, 7, 10};
@@ -60,8 +61,8 @@ void mvd()
    TGeoTranslation T;
    TGeoRotation R;
    TGeoCombiTrans *M;
-   char station_name[30], quadrant_name[30], carrier_name[30], heatsink_name[30], heatsinkpart_name[30];
-   int heatsinkpartno = 0;
+   char station_name[30], quadrant_name[30], carrier_name[30], heatsink_name[30], heatsinkpart_name[30], cable_name[30];
+   int heatsinkpartno, cablepartno = 0;
 
    // --- Sensor Assembly
    TGeoVolume *sensor = new TGeoVolumeAssembly("sensor");
@@ -112,33 +113,27 @@ void mvd()
       quadrant->AddNode(carrier, 1, M);
       // --- Heatsink
       R.SetAngles(0.,0.,0.);
-      // first part on carrier
+      // first part next to carrier
       sprintf(heatsinkpart_name, "heatsinkpart_%i", heatsinkpartno++);
       Float_t height_1 = heatsinkWidth[i]/2 - cd[1] + carrierClampOverlap - quadrantBeamOffset[i];
       Float_t width_1  = heatsinkWidth[i]/2 + quadrantBeamOffset[i];
-      Float_t thickness = (heatsinkThickness-cd[2])/2.;
+      Float_t thickness = heatsinkThickness;
       TGeoVolume* hs_part1 = manager->MakeBox(heatsinkpart_name, Al, width_1/2., height_1/2., thickness/2.);
       hs_part1->SetTransparency(10);
       hs_part1->SetLineColor(kGray);
-      T.SetTranslation(-width_1/2, -height_1/2-cd[1]+carrierClampOverlap, +(thickness+cd[2])/2*(1.0+explosion/0.2));
+      T.SetTranslation(-width_1/2, -height_1/2-cd[1]+carrierClampOverlap, 0);
       M = new TGeoCombiTrans(T,R);
       quadrant->AddNode(hs_part1, 1, M);
-      T.SetTranslation(-width_1/2, -height_1/2-cd[1]+carrierClampOverlap, -(thickness+cd[2])/2);
-      M = new TGeoCombiTrans(T,R);
-      quadrant->AddNode(hs_part1, 2, M);
-      // second part on carrier
+      // second part next to carrier
       sprintf(heatsinkpart_name, "heatsinkpart_%i", heatsinkpartno++);
       Float_t height_2 = heatsinkWidth[i] - height_1 - width_1;
       Float_t width_2  = heatsinkWidth[i]/2 - cd[0] + carrierClampOverlap + quadrantBeamOffset[i];
       TGeoVolume* hs_part2 = manager->MakeBox(heatsinkpart_name, Al, width_2/2, height_2/2, thickness/2);
       hs_part2->SetTransparency(10);
       hs_part2->SetLineColor(kGray);
-      T.SetTranslation(-width_2/2-cd[0]+carrierClampOverlap, -height_2/2, +(thickness+cd[2])/2);
+      T.SetTranslation(-width_2/2-cd[0]+carrierClampOverlap, -height_2/2, 0);
       M = new TGeoCombiTrans(T,R);
       quadrant->AddNode(hs_part2, 1, M);
-      T.SetTranslation(-width_2/2-cd[0]+carrierClampOverlap, -height_2/2, -(thickness+cd[2])/2);
-      M = new TGeoCombiTrans(T,R);
-      quadrant->AddNode(hs_part2, 2, M);
       // element to fill top and bottom if needed
       sprintf(heatsinkpart_name, "heatsinkpart_%i", heatsinkpartno++);
       Float_t height = (heatsinkHeight-heatsinkWidth[i])/2;
@@ -159,14 +154,22 @@ void mvd()
       // --- Sensors
       cout << "  # of sensors per quadrant: " << sensorCols[i]*sensorRows[i] << endl;
       R.SetAngles(0.,0.,0.);
+      TGeoVolume *cable;
       for (int k=0; k<sensorCols[i]; k++)
          {
+         if (k%2 == 0)
+            {
+               sprintf(cable_name, "cable_%i", cablepartno++);
+               cable = manager->MakeBox(cable_name, FPC, (cd[0] - sensorDimensionsActive[0]*k)/2, fpcWidth/2, fpcThickness/2);
+               cable->SetTransparency(0);
+               cable->SetLineColor(kSpring-1);
+            }
          for (int l=0; l<sensorRows[i]; l++)
             {
-            Float_t x_offset = -sensorDimensionsActive[0]*k;
             Float_t y_offset = - sensorPitch*l;
             Float_t z_offset = cd[2]/2;
             R.SetAngles(0.,0.,0.);
+            Float_t x_offset = - sensorDimensionsActive[0]*k;
             if (l%2 == 0) // front side
                z_offset = -z_offset;
             if (l%2 == 1) // back side
@@ -177,6 +180,29 @@ void mvd()
             T.SetTranslation(x_offset, y_offset, z_offset);
             M = new TGeoCombiTrans(T,R);
             quadrant->AddNode(sensor, k*sensorRows[i] + l, M);
+            if (k%2 == 0)
+               {
+               // add cable
+               R.SetAngles(0.,0.,0.);
+               Float_t fpc_x_offset = -(cd[0] - sensorDimensionsActive[0]*k)/2 - sensorDimensionsActive[0]*k;
+               Float_t fpc_y_offset = -fpcWidth/2 - sensorPitch*(l+1) - sensorDimensionsPassive[1];
+               Float_t fpc_z_offset = z_offset;
+               if (l%2 == 0) // front side
+                  {
+                  fpc_z_offset -= sensorDimensionsActive[2];
+                  fpc_z_offset -= glueThickness;
+                  fpc_z_offset -= fpcThickness/2;
+                  }
+               if (l%2 == 1) // back side
+                  {
+                  fpc_z_offset += sensorDimensionsActive[2];
+                  fpc_z_offset += glueThickness;
+                  fpc_z_offset += fpcThickness/2;
+                  }
+               T.SetTranslation(fpc_x_offset, fpc_y_offset, fpc_z_offset);
+               M = new TGeoCombiTrans(T,R);
+               quadrant->AddNode(cable, l+1, M);
+               }
             }
          }
       for (int j=0; j<4; j++)