Belle II Software development
GeoCDCCreator.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <cdc/geometry/GeoCDCCreator.h>
10#include <cdc/geometry/CDCGeometryPar.h>
11#include <cdc/geometry/CDCGeoControlPar.h>
12#include <cdc/simulation/CDCSimControlPar.h>
13#include <cdc/simulation/CDCSensitiveDetector.h>
14#include <simulation/background/BkgSensitiveDetector.h>
15
16#include <geometry/CreatorFactory.h>
17#include <geometry/Materials.h>
18
19#include <cmath>
20#include <boost/format.hpp>
21
22#include <G4Material.hh>
23#include <G4ProductionCuts.hh>
24#include <G4Box.hh>
25#include <G4Tubs.hh>
26#include <G4Torus.hh>
27#include <G4Trd.hh>
28#include <G4SubtractionSolid.hh>
29#include <G4Region.hh>
30#include <G4VSolid.hh>
31
32#include <G4Polycone.hh>
33#include <G4Cons.hh>
34#include <G4Colour.hh>
35#include <G4LogicalVolume.hh>
36#include <G4PVPlacement.hh>
37#include <G4Transform3D.hh>
38#include <G4VisAttributes.hh>
39#include <G4RotationMatrix.hh>
40#include <G4UserLimits.hh>
41#include <iostream>
42
43using namespace std;
44using namespace boost;
45
46namespace Belle2 {
52 using namespace geometry;
53
54 namespace CDC {
55
61
62 //-----------------------------------------------------------------
63 // Implementation
64 //-----------------------------------------------------------------
65
67 {
68 // Set job control params. before sensitivedetector and gometry construction
71
72 m_logicalCDC = 0;
73 m_physicalCDC = 0;
74 m_VisAttributes.clear();
75 m_VisAttributes.push_back(new G4VisAttributes(false)); // for "invisible"
76 m_userLimits.clear();
77 }
78
79
81 {
82 delete m_sensitive;
84 for (BkgSensitiveDetector* sensitiveDetector : m_BkgSensitiveRib4)
85 delete sensitiveDetector;
86 for (G4VisAttributes* visAttr : m_VisAttributes) delete visAttr;
87 m_VisAttributes.clear();
88 for (G4UserLimits* userLimits : m_userLimits) delete userLimits;
89 m_userLimits.clear();
90 }
91
92 void GeoCDCCreator::createGeometry(const CDCGeometry& geo, G4LogicalVolume& topVolume, geometry::GeometryTypes)
93 {
94
95 m_sensitive = new CDCSensitiveDetector("CDCSensitiveDetector", (2 * 24)* CLHEP::eV, 10 * CLHEP::MeV);
96
97 const G4double realTemperture = (273.15 + 23.) * CLHEP::kelvin;
98 G4Material* medHelium = geometry::Materials::get("CDCHeGas");
99 G4Material* medEthane = geometry::Materials::get("CDCEthaneGas");
100 G4Material* medAluminum = geometry::Materials::get("Al");
101 G4Material* medTungsten = geometry::Materials::get("W");
102 G4Material* medCFRP = geometry::Materials::get("CFRP");
103 G4Material* medNEMA_G10_Plate = geometry::Materials::get("NEMA_G10_Plate");
104 G4Material* medGlue = geometry::Materials::get("CDCGlue");
105 G4Material* medAir = geometry::Materials::get("Air");
106
107 G4double h2odensity = 1.000 * CLHEP::g / CLHEP::cm3;
108 G4double a = 1.01 * CLHEP::g / CLHEP::mole;
109 G4Element* elH = new G4Element("Hydrogen", "H", 1., a);
110 a = 16.00 * CLHEP::g / CLHEP::mole;
111 G4Element* elO = new G4Element("Oxygen", "O", 8., a);
112 G4Material* medH2O = new G4Material("Water", h2odensity, 2);
113 medH2O->AddElement(elH, 2);
114 medH2O->AddElement(elO, 1);
115 G4Material* medCopper = geometry::Materials::get("Cu");
116 G4Material* medHV = geometry::Materials::get("CDCHVCable");
117 //G4Material* medFiber = geometry::Materials::get("CDCOpticalFiber");
118 //G4Material* medCAT7 = geometry::Materials::get("CDCCAT7");
119 //G4Material* medTRG = geometry::Materials::get("CDCOpticalFiberTRG");
120
121 // Total cross section
122 const double rmax_innerWall = geo.getFiducialRmin();
123 const double rmin_outerWall = geo.getFiducialRmax();
124 const double diameter_senseWire = geo.getSenseDiameter();
125 const double diameter_fieldWire = geo.getFieldDiameter();
126 const double num_senseWire = static_cast<double>(geo.getNSenseWires());
127 const double num_fieldWire = static_cast<double>(geo.getNFieldWires());
128 double totalCS = M_PI * (rmin_outerWall * rmin_outerWall - rmax_innerWall * rmax_innerWall);
129
130 // Sense wire cross section
131 double senseCS = M_PI * (diameter_senseWire / 2) * (diameter_senseWire / 2) * num_senseWire;
132
133 // Field wire cross section
134 double fieldCS = M_PI * (diameter_fieldWire / 2) * (diameter_fieldWire / 2) * num_fieldWire;
135
136 // Density
137 const double denHelium = medHelium->GetDensity() / 2.0;
138 const double denEthane = medEthane->GetDensity() / 2.0;
139 const double denAluminum = medAluminum->GetDensity() * (fieldCS / totalCS);
140 const double denTungsten = medTungsten->GetDensity() * (senseCS / totalCS);
141 const double density = denHelium + denEthane + denAluminum + denTungsten;
142 G4Material* cdcMed = new G4Material("CDCGasWire", density, 4, kStateGas, realTemperture);
143 cdcMed->AddMaterial(medHelium, denHelium / density);
144 cdcMed->AddMaterial(medEthane, denEthane / density);
145 cdcMed->AddMaterial(medTungsten, denTungsten / density);
146 cdcMed->AddMaterial(medAluminum, denAluminum / density);
147
148 G4Material* cdcMedGas = cdcMed;
149
152 // std::cout << gcp.getMaterialDefinitionMode() << std::endl;
153
154 // if (cdcgp.getMaterialDefinitionMode() == 2) {
155 if (gcp.getMaterialDefinitionMode() == 2) {
156 const double density2 = denHelium + denEthane;
157 cdcMedGas = new G4Material("CDCRealGas", density2, 2, kStateGas, realTemperture);
158 cdcMedGas->AddMaterial(medHelium, denHelium / density2);
159 cdcMedGas->AddMaterial(medEthane, denEthane / density2);
160 }
161
162 if (gcp.getPrintMaterialTable()) {
163 G4cout << *(G4Material::GetMaterialTable());
164 }
165
166 const auto& mother = geo.getMotherVolume();
167 const auto& motherRmin = mother.getRmin();
168 const auto& motherRmax = mother.getRmax();
169 const auto& motherZ = mother.getZ();
170 G4Polycone* solid_cdc =
171 new G4Polycone("solidCDC", 0 * CLHEP::deg, 360.* CLHEP::deg,
172 mother.getNNodes(), motherZ.data(),
173 motherRmin.data(), motherRmax.data());
174 m_logicalCDC = new G4LogicalVolume(solid_cdc, medAir, "logicalCDC", 0, 0, 0);
175 m_physicalCDC = new G4PVPlacement(0, G4ThreeVector(geo.getGlobalOffsetX() * CLHEP::cm,
176 geo.getGlobalOffsetY() * CLHEP::cm,
177 geo.getGlobalOffsetZ() * CLHEP::cm), m_logicalCDC,
178 "physicalCDC", &topVolume, false, 0);
179
180 // Set up region for production cuts
181 G4Region* aRegion = new G4Region("CDCEnvelope");
182 m_logicalCDC->SetRegion(aRegion);
183 aRegion->AddRootLogicalVolume(m_logicalCDC);
184
185 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(0., 1., 0.)));
186 for (const auto& wall : geo.getOuterWalls()) {
187 const int iOuterWall = wall.getId();
188 const string wallName = wall.getName();
189 const double wallRmin = wall.getRmin();
190 const double wallRmax = wall.getRmax();
191 const double wallZfwd = wall.getZfwd();
192 const double wallZbwd = wall.getZbwd();
193 const double length = (wallZfwd - wallZbwd) / 2.0;
194
195
196 G4Material* medWall;
197 if (strstr((wallName).c_str(), "MiddleWall") != nullptr) {
198 medWall = medCFRP;
199 } else {
200 medWall = medAluminum;
201 }
202 G4Tubs* outerWallTubeShape = new G4Tubs("solid" + wallName, wallRmin * CLHEP::cm,
203 wallRmax * CLHEP::cm, length * CLHEP::cm, 0 * CLHEP::deg, 360.*CLHEP::deg);
204
205 G4LogicalVolume* outerWallTube = new G4LogicalVolume(outerWallTubeShape, medWall, "solid" + wallName, 0, 0, 0);
206 outerWallTube->SetVisAttributes(m_VisAttributes.back());
207 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (length + wallZbwd)*CLHEP::cm), outerWallTube, "logical" + wallName,
208 m_logicalCDC, false, iOuterWall);
209 }
210
211
212 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(0., 1., 0.)));
213 for (const auto& wall : geo.getInnerWalls()) {
214 const string wallName = wall.getName();
215 const double wallRmin = wall.getRmin();
216 const double wallRmax = wall.getRmax();
217 const double wallZfwd = wall.getZfwd();
218 const double wallZbwd = wall.getZbwd();
219 const double length = (wallZfwd - wallZbwd) / 2.0;
220 const int iInnerWall = wall.getId();
221
222 G4Material* medWall;
223 if (strstr(wallName.c_str(), "MiddleWall") != nullptr) {
224 medWall = medCFRP;
225 } else if (strstr(wallName.c_str(), "MiddleGlue") != nullptr) { // Glue layer 0.005 mmt
226 medWall = medGlue;
227 } else { // Al layer 0.1 mmt
228 medWall = medAluminum;
229 }
230
231 G4Tubs* innerWallTubeShape = new G4Tubs("solid" + wallName, wallRmin * CLHEP::cm,
232 wallRmax * CLHEP::cm, length * CLHEP::cm, 0 * CLHEP::deg, 360.*CLHEP::deg);
233 G4LogicalVolume* innerWallTube = new G4LogicalVolume(innerWallTubeShape, medWall, "logical" + wallName, 0, 0, 0);
234 innerWallTube->SetVisAttributes(m_VisAttributes.back());
235 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (length + wallZbwd)*CLHEP::cm), innerWallTube, "physical" + wallName,
236 m_logicalCDC, false, iInnerWall);
237
238
239 }
240
241
242
243 //
244 // Construct sensitive layers.
245 //
246 const int nSLayer = geo.getNSenseLayers();
247 const double length_feedthrough = geo.getFeedthroughLength();
248 for (int iSLayer = 0; iSLayer < nSLayer; ++iSLayer) {
249 const auto& endplate = geo.getEndPlate(iSLayer);
250 const int nEPLayer = endplate.getNEndPlateLayers();
251 // Get parameters for sensitive layer: left, middle and right.
252 double rmin_sensitive_left, rmax_sensitive_left;
253 double rmin_sensitive_middle, rmax_sensitive_middle;
254 double rmin_sensitive_right, rmax_sensitive_right;
255 double zback_sensitive_left, zfor_sensitive_left;
256 double zback_sensitive_middle, zfor_sensitive_middle;
257 double zback_sensitive_right, zfor_sensitive_right;
258
259 if (iSLayer == 0) {
260 const auto& epLayerBwd = endplate.getEndPlateLayer(0);
261 const auto& epLayerFwd = endplate.getEndPlateLayer(nEPLayer / 2);
262 const auto& senseLayer = geo.getSenseLayer(iSLayer);
263 const auto& fieldLayer = geo.getFieldLayer(iSLayer);
264
265 rmin_sensitive_left = epLayerBwd.getRmax();
266 rmax_sensitive_left = fieldLayer.getR();
267 zback_sensitive_left = senseLayer.getZbwd();
268 zfor_sensitive_left = epLayerBwd.getZfwd();
269
270 rmin_sensitive_middle = (geo.getInnerWall(0)).getRmax();
271 rmax_sensitive_middle = fieldLayer.getR();
272 zback_sensitive_middle = epLayerBwd.getZfwd();
273 zfor_sensitive_middle = epLayerFwd.getZbwd();
274
275 rmin_sensitive_right = epLayerFwd.getRmax();
276 rmax_sensitive_right = fieldLayer.getR();
277 zback_sensitive_right = epLayerFwd.getZbwd();
278 zfor_sensitive_right = senseLayer.getZfwd();
279 } else if (iSLayer >= 1 && iSLayer <= 14) {
280 const auto& epLayerBwd = endplate.getEndPlateLayer(1);
281 const auto& epLayerFwd = endplate.getEndPlateLayer((nEPLayer / 2) + 1);
282 const auto& senseLayer = geo.getSenseLayer(iSLayer);
283 const auto& fieldLayerIn = geo.getFieldLayer(iSLayer - 1);
284 const auto& fieldLayerOut = geo.getFieldLayer(iSLayer);
285
286 rmin_sensitive_left = epLayerBwd.getRmax();
287 rmax_sensitive_left = fieldLayerOut.getR();
288 zback_sensitive_left = senseLayer.getZbwd();
289 zfor_sensitive_left = epLayerBwd.getZfwd();
290
291 rmin_sensitive_middle = fieldLayerIn.getR();
292 rmax_sensitive_middle = fieldLayerOut.getR();
293 zback_sensitive_middle = epLayerBwd.getZfwd();
294 zfor_sensitive_middle = epLayerFwd.getZbwd();
295
296 rmin_sensitive_right = epLayerFwd.getRmax();
297 rmax_sensitive_right = fieldLayerOut.getR();
298 zback_sensitive_right = epLayerFwd.getZbwd();
299 zfor_sensitive_right = senseLayer.getZfwd();
300 } else if (iSLayer >= 15 && iSLayer <= 18) {
301 const auto& epLayerBwd = endplate.getEndPlateLayer(1);
302 const auto& epLayerFwd = endplate.getEndPlateLayer(nEPLayer / 2);
303 const auto& senseLayer = geo.getSenseLayer(iSLayer);
304 const auto& fieldLayerIn = geo.getFieldLayer(iSLayer - 1);
305 const auto& fieldLayerOut = geo.getFieldLayer(iSLayer);
306
307 rmin_sensitive_left = epLayerBwd.getRmax();
308 rmax_sensitive_left = fieldLayerOut.getR();
309 zback_sensitive_left = senseLayer.getZbwd();
310 zfor_sensitive_left = epLayerBwd.getZfwd();
311
312 rmin_sensitive_middle = fieldLayerIn.getR();
313 rmax_sensitive_middle = fieldLayerOut.getR();
314 zback_sensitive_middle = epLayerBwd.getZfwd();
315 zfor_sensitive_middle = epLayerFwd.getZbwd();
316
317 rmin_sensitive_right = epLayerFwd.getRmax();
318 rmax_sensitive_right = fieldLayerOut.getR();
319 zback_sensitive_right = epLayerFwd.getZbwd();
320 zfor_sensitive_right = senseLayer.getZfwd();
321 } else if (iSLayer >= 19 && iSLayer < 55) {
322 const auto& epLayerBwd = endplate.getEndPlateLayer(0);
323 const auto& epLayerFwd = endplate.getEndPlateLayer(nEPLayer / 2);
324 const auto& senseLayer = geo.getSenseLayer(iSLayer);
325 const auto& fieldLayerIn = geo.getFieldLayer(iSLayer - 1);
326 const auto& fieldLayerOut = geo.getFieldLayer(iSLayer);
327
328 rmin_sensitive_left = epLayerBwd.getRmax();
329 rmax_sensitive_left = fieldLayerOut.getR();
330 zback_sensitive_left = senseLayer.getZbwd();
331 zfor_sensitive_left = epLayerBwd.getZfwd();
332
333 rmin_sensitive_middle = fieldLayerIn.getR();
334 rmax_sensitive_middle = fieldLayerOut.getR();
335 zback_sensitive_middle = epLayerBwd.getZfwd();
336 zfor_sensitive_middle = epLayerFwd.getZbwd();
337
338 rmin_sensitive_right = epLayerFwd.getRmax();
339 rmax_sensitive_right = fieldLayerOut.getR();
340 zback_sensitive_right = epLayerFwd.getZbwd();
341 zfor_sensitive_right = senseLayer.getZfwd();
342
343 } else if (iSLayer == 55) {
344
345 const auto& epLayerBwdIn = endplate.getEndPlateLayer(0);
346 const auto& epLayerBwdOut = endplate.getEndPlateLayer((nEPLayer / 2) - 1);
347 const auto& epLayerFwdIn = endplate.getEndPlateLayer(nEPLayer / 2);
348 const auto& epLayerFwdOut = endplate.getEndPlateLayer(nEPLayer - 1);
349 const auto& senseLayer = geo.getSenseLayer(iSLayer);
350 // const auto& fieldLayerIn = geo.getFieldLayer(iSLayer - 1);
351 int iSLayerMinus1 = iSLayer - 1; //avoid cpp-check warning
352 const auto& fieldLayerIn = geo.getFieldLayer(iSLayerMinus1); //avoid cpp-check warning
353 rmin_sensitive_left = epLayerBwdIn.getRmax();
354 rmax_sensitive_left = epLayerBwdOut.getRmax();
355 zback_sensitive_left = senseLayer.getZbwd();
356 zfor_sensitive_left = epLayerBwdIn.getZfwd();
357
358 rmin_sensitive_middle = fieldLayerIn.getR();
359 rmax_sensitive_middle = (geo.getOuterWall(0)).getRmin();
360 zback_sensitive_middle = epLayerBwdIn.getZfwd();
361 zfor_sensitive_middle = epLayerFwdIn.getZbwd();
362
363 rmin_sensitive_right = epLayerFwdIn.getRmax();
364 rmax_sensitive_right = epLayerFwdOut.getRmax();
365 zback_sensitive_right = epLayerFwdIn.getZbwd();
366 zfor_sensitive_right = senseLayer.getZfwd();
367
368 } else {
369 B2ERROR("Undefined sensitive layer : " << iSLayer);
370 continue;
371 }
372
373
374 // Check if build left sensitive tube
375 if ((zfor_sensitive_left - zback_sensitive_left) > length_feedthrough) {
376 // std::cout <<"left doif " << iSLayer <<" "<< zfor_sensitive_left - zback_sensitive_left << std::endl;
377 //==========================================================
378 // zback_sensitive_left
379 // |
380 // \|/
381 // _____________________
382 // | |// 1 // | |
383 // | |==ft====|2 | (ft = feedthrouth)
384 // |_______|____1___|__|
385 // |_______|___________|
386 // |
387 // \|/
388 // zfor_sensitive_left
389 //==========================================================
390
391 // Build a tube with metarial cdcMed for area 1
392 G4Tubs* leftTubeShape = new G4Tubs((format("solidCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), rmin_sensitive_left * CLHEP::cm,
393 rmax_sensitive_left * CLHEP::cm, length_feedthrough * CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
394 G4LogicalVolume* leftTube = new G4LogicalVolume(leftTubeShape, cdcMed,
395 (format("logicalCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), 0, 0, 0);
396 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zback_sensitive_left + length_feedthrough / 2.0)*CLHEP::cm), leftTube,
397 (format("physicalCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
398 // Build left sensitive tube (area 2)
399 G4Tubs* leftSensitiveTubeShape = new G4Tubs((format("solidSD_CDCLayer_%1%_left") % iSLayer).str().c_str(),
400 rmin_sensitive_left * CLHEP::cm, rmax_sensitive_left * CLHEP::cm,
401 (zfor_sensitive_left - zback_sensitive_left - length_feedthrough)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
402 G4LogicalVolume* leftSensitiveTube = new G4LogicalVolume(leftSensitiveTubeShape, cdcMed,
403 (format("logicalSD_CDCLayer_%1%_left") % iSLayer).str().c_str(), 0, 0, 0);
404 leftSensitiveTube->SetSensitiveDetector(m_sensitive);
405 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_left + zback_sensitive_left + length_feedthrough)*CLHEP::cm / 2.0),
406 leftSensitiveTube, (format("physicalSD_CDCLayer_%1%_left") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
407 } else {
408 // std::cout <<"left doelse " << iSLayer << std::endl;
409 //==========================================================
410 // zback_sensitive_left
411 // |
412 // \|/
413 // _________________________
414 // | |//// 1 ////| 2 |
415 // | |======ft======== (ft = feedthrouth)
416 // |_______|____1______| 2 |
417 // |_______|___________|___|
418 // |
419 // \|/
420 // zfor_sensitive_left
421 //==========================================================
422
423 // Build a tube with metarial cdcMed for area 1
424 G4Tubs* leftTubeShape = new G4Tubs((format("solidCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), rmin_sensitive_left * CLHEP::cm,
425 rmax_sensitive_left * CLHEP::cm, (zfor_sensitive_left - zback_sensitive_left)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
426 G4LogicalVolume* leftTube = new G4LogicalVolume(leftTubeShape, cdcMed,
427 (format("logicalCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), 0, 0, 0);
428 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_left + zback_sensitive_left)*CLHEP::cm / 2.0), leftTube,
429 (format("physicalCDCLayer_%1%_leftTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
430
431
432 // Build a tube with metarial cdcMed for area 2
433 G4Tubs* leftMidTubeShape = new G4Tubs((format("solidCDCLayer_%1%_leftMidTube") % iSLayer).str().c_str(),
434 rmin_sensitive_middle * CLHEP::cm, rmax_sensitive_middle * CLHEP::cm,
435 (length_feedthrough - zfor_sensitive_left + zback_sensitive_left)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
436 G4LogicalVolume* leftMidTube = new G4LogicalVolume(leftMidTubeShape, cdcMed,
437 (format("logicalCDCLayer_%1%_leftMidTube") % iSLayer).str().c_str(), 0, 0, 0);
438
439 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (length_feedthrough + zfor_sensitive_left + zback_sensitive_left)*CLHEP::cm / 2.0),
440 leftMidTube, (format("physicalCDCLayer_%1%_leftMidTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
441
442 // Reset zback_sensitive_middle
443 zback_sensitive_middle = length_feedthrough + zback_sensitive_left;
444 }
445
446 // Check if build right sensitive tube
447 if ((zfor_sensitive_right - zback_sensitive_right) > length_feedthrough) {
448 // std::cout <<"right doif" << iSLayer <<" "<< zfor_sensitive_right - zback_sensitive_right << std::endl;
449 //==========================================================
450 // zfor_sensitive_right
451 // |
452 // \|/
453 // _____________________________
454 // | | 1 |///////|
455 // | 2 |====ft=====|///////| (ft = feedthrouth)
456 // |_______|____1______|_______|
457 // |_______|___________|_______|
458 // |
459 // \|/
460 // zback_sensitive_right
461 //==========================================================
462
463 // Build a tube with metarial cdcMed for area 1
464 G4Tubs* rightTubeShape = new G4Tubs((format("solidCDCLayer_%1%_rightTube") % iSLayer).str().c_str(),
465 rmin_sensitive_right * CLHEP::cm, rmax_sensitive_right * CLHEP::cm, length_feedthrough * CLHEP::cm / 2.0, 0 * CLHEP::deg,
466 360.*CLHEP::deg);
467 G4LogicalVolume* rightTube = new G4LogicalVolume(rightTubeShape, cdcMed,
468 (format("logicalCDCLayer_%1%_rightTube") % iSLayer).str().c_str(), 0, 0, 0);
469
470 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_right - length_feedthrough / 2.0)*CLHEP::cm), rightTube,
471 (format("physicalCDCLayer_%1%_rightTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
472
473
474 // Build right sensitive tube (area 2)
475 G4Tubs* rightSensitiveTubeShape = new G4Tubs((format("solidSD_CDCLayer_%1%_right") % iSLayer).str().c_str(),
476 rmin_sensitive_right * CLHEP::cm, rmax_sensitive_right * CLHEP::cm,
477 (zfor_sensitive_right - zback_sensitive_right - length_feedthrough)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
478 G4LogicalVolume* rightSensitiveTube = new G4LogicalVolume(rightSensitiveTubeShape, cdcMed,
479 (format("logicalSD_CDCLayer_%1%_right") % iSLayer).str().c_str(), 0, 0, 0);
480 rightSensitiveTube->SetSensitiveDetector(m_sensitive);
481
482 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_right + zback_sensitive_right - length_feedthrough)*CLHEP::cm / 2.0),
483 rightSensitiveTube, (format("physicalSD_CDCLayer_%1%_right") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
484
485 } else {
486 // std::cout <<"right doelse" << iSLayer << std::endl;
487 //==========================================================
488 // zfor_sensitive_right
489 // |
490 // \|/
491 // _____________________________
492 // | | 1 |///////|
493 // |============ft=====|///////| (ft = feedthrouth)
494 // | |____1______|_______|
495 // |_______|___________|_______|
496 // |
497 // \|/
498 // zback_sensitive_right
499 //==========================================================
500
501 // Build a tube with metarial cdcMed for area 1
502 G4Tubs* rightTubeShape = new G4Tubs((format("solidCDCLayer_%1%_rightTube") % iSLayer).str().c_str(),
503 rmin_sensitive_right * CLHEP::cm, rmax_sensitive_right * CLHEP::cm, (zfor_sensitive_right - zback_sensitive_right)*CLHEP::cm / 2.0,
504 0 * CLHEP::deg, 360.*CLHEP::deg);
505 G4LogicalVolume* rightTube = new G4LogicalVolume(rightTubeShape, cdcMed,
506 (format("logicalCDCLayer_%1%_rightTube") % iSLayer).str().c_str(), 0, 0, 0);
507
508 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_right + zback_sensitive_right)*CLHEP::cm / 2.0), rightTube,
509 (format("physicalCDCLayer_%1%_rightTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
510
511
512 // Build a tube with metarial cdcMed for area 2
513 G4Tubs* rightMidTubeShape = new G4Tubs((format("solidCDCLayer_%1%_rightMidTube") % iSLayer).str().c_str(),
514 rmin_sensitive_middle * CLHEP::cm, rmax_sensitive_middle * CLHEP::cm,
515 (length_feedthrough - zfor_sensitive_right + zback_sensitive_right)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
516 G4LogicalVolume* rightMidTube = new G4LogicalVolume(rightMidTubeShape, cdcMed,
517 (format("logicalCDCLayer_%1%_rightMidTube") % iSLayer).str().c_str(), 0, 0, 0);
518 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zback_sensitive_right - length_feedthrough + zfor_sensitive_right)*CLHEP::cm / 2.0),
519 rightMidTube, (format("physicalCDCLayer_%1%_rightMidTube") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
520
521 // Reset zback_sensitive_middle
522 zfor_sensitive_middle = zfor_sensitive_right - length_feedthrough;
523 }
524
525
526 // Middle sensitive tube
527 G4Tubs* middleSensitiveTubeShape = new G4Tubs((format("solidSD_CDCLayer_%1%_middle") % iSLayer).str().c_str(),
528 rmin_sensitive_middle * CLHEP::cm, rmax_sensitive_middle * CLHEP::cm,
529 (zfor_sensitive_middle - zback_sensitive_middle)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
530 G4LogicalVolume* middleSensitiveTube = new G4LogicalVolume(middleSensitiveTubeShape, cdcMedGas,
531 (format("logicalSD_CDCLayer_%1%_middle") % iSLayer).str().c_str(), 0, 0, 0);
532 //hard-coded temporarily
533 //need to create an object per layer ??? to be checked later
534 G4UserLimits* uLimits = new G4UserLimits(8.5 * CLHEP::cm);
535 m_userLimits.push_back(uLimits);
536 middleSensitiveTube->SetUserLimits(uLimits);
537 middleSensitiveTube->SetSensitiveDetector(m_sensitive);
538
539 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfor_sensitive_middle + zback_sensitive_middle)*CLHEP::cm / 2.0), middleSensitiveTube,
540 (format("physicalSD_CDCLayer_%1%_middle") % iSLayer).str().c_str(), m_logicalCDC, false, iSLayer);
541
542 // if (cdcgp.getMaterialDefinitionMode() == 2) {
543 if (gcp.getMaterialDefinitionMode() == 2) {
544 G4String sName = "sWire";
545 const G4int jc = 0;
546 B2Vector3D wb0 = cdcgp.wireBackwardPosition(iSLayer, jc);
547 // G4double rsense0 = wb0.Perp();
548 B2Vector3D wf0 = cdcgp.wireForwardPosition(iSLayer, jc);
549 G4double tAtZ0 = -wb0.Z() / (wf0.Z() - wb0.Z()); //t: param. along the wire
550 B2Vector3D wAtZ0 = wb0 + tAtZ0 * (wf0 - wb0);
551 //additional chop of wire; must be larger than 126um*(1/10), where 126um is the field wire diameter; 1/10: approx. stereo angle
552 const G4double epsl = 25.e-4; // (in cm);
553 G4double reductionBwd = (zback_sensitive_middle + epsl) / wb0.Z();
554 //chop the wire at zback_sensitive_middle for avoiding overlap; this is because the wire length defined by wb0 and wf0 is larger than the length of middle sensitive tube
555 wb0 = reductionBwd * (wb0 - wAtZ0) + wAtZ0;
556 //chop wire at zfor_sensitive_middle for avoiding overlap
557 G4double reductionFwd = (zfor_sensitive_middle - epsl) / wf0.Z();
558 wf0 = reductionFwd * (wf0 - wAtZ0) + wAtZ0;
559
560 const G4double wireHalfLength = 0.5 * (wf0 - wb0).Mag() * CLHEP::cm;
561 const G4double sWireRadius = 0.5 * cdcgp.senseWireDiameter() * CLHEP::cm;
562 // const G4double sWireRadius = 15.e-4 * CLHEP::cm;
563 G4Tubs* middleSensitiveSwireShape = new G4Tubs(sName, 0., sWireRadius, wireHalfLength, 0., 360. * CLHEP::deg);
564 G4LogicalVolume* middleSensitiveSwire = new G4LogicalVolume(middleSensitiveSwireShape, medTungsten, sName);
565 // middleSensitiveSwire->SetSensitiveDetector(m_sensitive);
566 middleSensitiveSwire->SetVisAttributes(m_VisAttributes.front()); // <- to speed up visualization
567
568 G4String fName = "fWire";
569 const G4double fWireRadius = 0.5 * cdcgp.fieldWireDiameter() * CLHEP::cm;
570 G4Tubs* middleSensitiveFwireShape = new G4Tubs(fName, 0., fWireRadius, wireHalfLength, 0., 360. * CLHEP::deg);
571 G4LogicalVolume* middleSensitiveFwire = new G4LogicalVolume(middleSensitiveFwireShape, medAluminum, fName);
572 // middleSensitiveFwire->SetSensitiveDetector(m_sensitive);
573 middleSensitiveFwire->SetVisAttributes(m_VisAttributes.front()); // <- to speed up visualization
574
575 const G4double diameter = cdcgp.fieldWireDiameter();
576
577 const G4int nCells = cdcgp.nWiresInLayer(iSLayer);
578 const G4double dphi = M_PI / nCells;
579 const B2Vector3D unitZ(0., 0., 1.);
580
581 for (int ic = 0; ic < nCells; ++ic) {
582 //define sense wire
583 B2Vector3D wb = cdcgp.wireBackwardPosition(iSLayer, ic);
584 B2Vector3D wf = cdcgp.wireForwardPosition(iSLayer, ic);
585 G4double tAtZ02 = -wb.Z() / (wf.Z() - wb.Z());
586 B2Vector3D wAtZ02 = wb + tAtZ02 * (wf - wb);
587 G4double reductionBwd2 = (zback_sensitive_middle + epsl) / wb.Z();
588 wb = reductionBwd2 * (wb - wAtZ02) + wAtZ02;
589 G4double reductionFwd2 = (zfor_sensitive_middle - epsl) / wf.Z();
590 wf = reductionFwd2 * (wf - wAtZ02) + wAtZ02;
591
592 G4double thetaYZ = -asin((wf - wb).Y() / (wf - wb).Mag());
593
594 B2Vector3D fMinusBInZX((wf - wb).X(), 0., (wf - wb).Z());
595 G4double thetaZX = asin((unitZ.Cross(fMinusBInZX)).Y() / fMinusBInZX.Mag());
596 G4RotationMatrix rotM;
597 // std::cout <<"deg,rad= " << CLHEP::deg <<" "<< CLHEP::rad << std::endl;
598 rotM.rotateX(thetaYZ * CLHEP::rad);
599 rotM.rotateY(thetaZX * CLHEP::rad);
600
601 G4ThreeVector xyz(0.5 * (wb.X() + wf.X()) * CLHEP::cm,
602 0.5 * (wb.Y() + wf.Y()) * CLHEP::cm, 0.);
603
604 // std::cout <<"0 x,y= " << xyz.getX() <<" "<< xyz.getY() << std::endl;
605 //Calling G4PVPlacement with G4Transform3D is convenient because it rotates the object instead of rotating the coordinate-frame; rotM is copied so it does not have to be created on heep by new.
606 new G4PVPlacement(G4Transform3D(rotM, xyz), middleSensitiveSwire, sName, middleSensitiveTube, false, ic);
607
608 //define field wire #1 (placed at the same phi but at the outer r boundary)
609 B2Vector3D wbF = wb;
610 G4double rF = rmax_sensitive_middle - 0.5 * diameter;
611 // std::cout <<"iSLayer,rF= " << iSLayer <<" "<< rF <<" "<< std::endl;
612 G4double phi = atan2(wbF.Y(), wbF.X());
613 wbF.SetX(rF * cos(phi));
614 wbF.SetY(rF * sin(phi));
615
616 B2Vector3D wfF = wf;
617 rF = rmax_sensitive_middle - 0.5 * diameter;
618 phi = atan2(wfF.Y(), wfF.X());
619 wfF.SetX(rF * cos(phi));
620 wfF.SetY(rF * sin(phi));
621
622 thetaYZ = -asin((wfF - wbF).Y() / (wfF - wbF).Mag());
623
624 fMinusBInZX = wfF - wbF;
625 fMinusBInZX.SetY(0.);
626 thetaZX = asin((unitZ.Cross(fMinusBInZX)).Y() / fMinusBInZX.Mag());
627
628 G4RotationMatrix rotM1;
629 rotM1.rotateX(thetaYZ * CLHEP::rad);
630 rotM1.rotateY(thetaZX * CLHEP::rad);
631
632 xyz.setX(0.5 * (wbF.X() + wfF.X()) * CLHEP::cm);
633 xyz.setY(0.5 * (wbF.Y() + wfF.Y()) * CLHEP::cm);
634
635 if (iSLayer != nSLayer - 1) {
636 // std::cout <<"1 x,y= " << xyz.getX() <<" "<< xyz.getY() << std::endl;
637 new G4PVPlacement(G4Transform3D(rotM1, xyz), middleSensitiveFwire, fName, middleSensitiveTube, false, ic);
638 }
639
640 //define field wire #2 (placed at the same radius but shifted by dphi)
641 wbF = wb;
642 rF = wbF.Perp();
643 phi = atan2(wbF.Y(), wbF.X());
644 wbF.SetX(rF * cos(phi + dphi));
645 wbF.SetY(rF * sin(phi + dphi));
646
647 wfF = wf;
648 rF = wfF.Perp();
649 phi = atan2(wfF.Y(), wfF.X());
650 wfF.SetX(rF * cos(phi + dphi));
651 wfF.SetY(rF * sin(phi + dphi));
652
653 thetaYZ = -asin((wfF - wbF).Y() / (wfF - wbF).Mag());
654
655 fMinusBInZX = wfF - wbF;
656 fMinusBInZX.SetY(0.);
657 thetaZX = asin((unitZ.Cross(fMinusBInZX)).Y() / fMinusBInZX.Mag());
658
659 G4RotationMatrix rotM2;
660 rotM2.rotateX(thetaYZ * CLHEP::rad);
661 rotM2.rotateY(thetaZX * CLHEP::rad);
662
663 xyz.setX(0.5 * (wbF.X() + wfF.X()) * CLHEP::cm);
664 xyz.setY(0.5 * (wbF.Y() + wfF.Y()) * CLHEP::cm);
665
666 // std::cout <<"2 x,y= " << xyz.getX() <<" "<< xyz.getY() << std::endl;
667 new G4PVPlacement(G4Transform3D(rotM2, xyz), middleSensitiveFwire, fName, middleSensitiveTube, false, ic + nCells);
668
669 //define field wire #3 (placed at the cell corner)
670 wbF = wb;
671 rF = rmax_sensitive_middle - 0.5 * diameter;
672 phi = atan2(wbF.Y(), wbF.X());
673 wbF.SetX(rF * cos(phi + dphi));
674 wbF.SetY(rF * sin(phi + dphi));
675
676 wfF = wf;
677 rF = rmax_sensitive_middle - 0.5 * diameter;
678 phi = atan2(wfF.Y(), wfF.X());
679 wfF.SetX(rF * cos(phi + dphi));
680 wfF.SetY(rF * sin(phi + dphi));
681
682 thetaYZ = -asin((wfF - wbF).Y() / (wfF - wbF).Mag());
683
684 fMinusBInZX = wfF - wbF;
685 fMinusBInZX.SetY(0.);
686 thetaZX = asin((unitZ.Cross(fMinusBInZX)).Y() / fMinusBInZX.Mag());
687
688 G4RotationMatrix rotM3;
689 rotM3.rotateX(thetaYZ * CLHEP::rad);
690 rotM3.rotateY(thetaZX * CLHEP::rad);
691
692 xyz.setX(0.5 * (wbF.X() + wfF.X()) * CLHEP::cm);
693 xyz.setY(0.5 * (wbF.Y() + wfF.Y()) * CLHEP::cm);
694
695 if (iSLayer != nSLayer - 1) {
696 new G4PVPlacement(G4Transform3D(rotM3, xyz), middleSensitiveFwire, fName, middleSensitiveTube, false, ic + 2 * nCells);
697 }
698 } // end of wire loop
699 } // end of wire definitions
700
701 }
702 //
703 // Endplates.
704 //
705
706 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(1., 1., 0.)));
707 for (const auto& endplate : geo.getEndPlates()) {
708 for (const auto& epLayer : endplate.getEndPlateLayers()) {
709 const int iEPLayer = epLayer.getILayer();
710 const string name = epLayer.getName();
711 const double rmin = epLayer.getRmin();
712 const double rmax = epLayer.getRmax();
713 const double zbwd = epLayer.getZbwd();
714 const double zfwd = epLayer.getZfwd();
715 const double length = (zfwd - zbwd) / 2.0;
716
717 G4Tubs* tube = new G4Tubs("solidCDCEndplate" + name, rmin * CLHEP::cm,
718 rmax * CLHEP::cm, length * CLHEP::cm, 0 * CLHEP::deg, 360.*CLHEP::deg);
719 G4LogicalVolume* logical = new G4LogicalVolume(tube, Materials::get("G4_Al"),
720 "logicalCDCEndplate" + name, 0, 0);
721 logical->SetVisAttributes(m_VisAttributes.back());
722 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (zfwd + zbwd)*CLHEP::cm / 2.0), logical,
723 "physicalCDCEndplate" + name, m_logicalCDC, false, iEPLayer);
724
725 }
726 }
727
728
729 // Construct electronics boards
730 for (const auto& frontend : geo.getFrontends()) {
731
732 const int iEB = frontend.getId();
733 const double ebInnerR = frontend.getRmin();
734 const double ebOuterR = frontend.getRmax();
735 const double ebBZ = frontend.getZbwd();
736 const double ebFZ = frontend.getZfwd();
737
738 G4Tubs* ebTubeShape = new G4Tubs((format("solidSD_ElectronicsBoard_Layer%1%") % iEB).str().c_str(), ebInnerR * CLHEP::cm,
739 ebOuterR * CLHEP::cm, (ebFZ - ebBZ)*CLHEP::cm / 2.0, 0 * CLHEP::deg, 360.*CLHEP::deg);
740
741 G4LogicalVolume* ebTube = new G4LogicalVolume(ebTubeShape, medNEMA_G10_Plate,
742 (format("logicalSD_ElectronicsBoard_Layer%1%") % iEB).str().c_str(), 0, 0, 0);
743 if (!m_bkgsensitive) m_bkgsensitive = new BkgSensitiveDetector("CDC", iEB);
744 ebTube->SetSensitiveDetector(m_bkgsensitive);
745 ebTube->SetVisAttributes(m_VisAttributes.back());
746 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (ebFZ + ebBZ)*CLHEP::cm / 2.0), ebTube,
747 (format("physicalSD_ElectronicsBoard_Layer%1%") % iEB).str().c_str(), m_logicalCDC, false, iEB);
748 }
749
750 //
751 // Construct neutron shield.
752 //
754
755 //
756 // construct covers.
757 //
758 createCovers(geo);
759
760 //
761 // construct covers.
762 //
763 createCover2s(geo);
764
765 //
766 // Construct ribs.
767 //
768 for (const auto& rib : geo.getRibs()) {
769
770 const int id = rib.getId();
771 const double length = rib.getLength();
772 const double width = rib.getWidth();
773 const double thick = rib.getThick();
774 const double rotx = rib.getRotX();
775 const double roty = rib.getRotY();
776 const double rotz = rib.getRotZ();
777 const double x = rib.getX();
778 const double y = rib.getY();
779 const double z = rib.getZ();
780 const int offset = rib.getOffset();
781 const int ndiv = rib.getNDiv();
782
783 const string solidName = "solidRib" + to_string(id);
784 const string logicalName = "logicalRib" + to_string(id);
785 G4Box* boxShape = new G4Box(solidName, 0.5 * length * CLHEP::cm,
786 0.5 * width * CLHEP::cm,
787 0.5 * thick * CLHEP::cm);
788
789 const double rmax = 0.5 * length;
790 const double rmin = max((rmax - thick), 0.);
791 G4Tubs* tubeShape = new G4Tubs(solidName,
792 rmin * CLHEP::cm,
793 rmax * CLHEP::cm,
794 0.5 * width * CLHEP::cm,
795 0.,
796 360. * CLHEP::deg);
797
798 //G4LogicalVolume* logicalV = new G4LogicalVolume(boxShape, medAluminum,
799 // logicalName, 0, 0, 0);
800 // ID depndent material definition, Aluminum is default
801 G4LogicalVolume* logicalV = new G4LogicalVolume(boxShape, medAluminum, logicalName, 0, 0, 0);
802 if (id > 39 && id < 78) // Cu
803 logicalV = new G4LogicalVolume(boxShape, medCopper, logicalName, 0, 0, 0);
804 if ((id > 77 && id < 94) || (id > 131 && id < 146)) // G10
805 logicalV = new G4LogicalVolume(boxShape, medNEMA_G10_Plate, logicalName, 0, 0, 0);
806 if (id > 93 && id < 110) // Cu
807 logicalV = new G4LogicalVolume(tubeShape, medCopper, logicalName, 0, 0, 0);
808 if (id > 109 && id < 126) // H2O
809 logicalV = new G4LogicalVolume(tubeShape, medH2O, logicalName, 0, 0, 0);
810 if (id > 127 && id < 132) // HV
811 logicalV = new G4LogicalVolume(boxShape, medHV, logicalName, 0, 0, 0);
812 /*if( id > 145 && id < 149 )// Fiber
813 logicalV = new G4LogicalVolume(boxShape, medFiber, logicalName, 0, 0, 0);
814 if( id > 148 && id < 158 )// Fiber
815 logicalV = new G4LogicalVolume(boxShape, medCAT7, logicalName, 0, 0, 0);
816 if( id > 157 && id < 164 )// Fiber
817 logicalV = new G4LogicalVolume(boxShape, medTRG, logicalName, 0, 0, 0);*/
818
819 logicalV->SetVisAttributes(m_VisAttributes.back());
820
821 const double phi = 360.0 / ndiv;
822
823 G4RotationMatrix rot = G4RotationMatrix();
824 double dz = thick;
825 if (id > 93 && id < 126) dz = 0;
826
827 G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - dz * CLHEP::cm / 2.0);
828 rot.rotateX(rotx);
829 rot.rotateY(roty);
830 rot.rotateZ(rotz);
831 if (offset) {
832 rot.rotateZ(0.5 * phi * CLHEP::deg);
833 arm.rotateZ(0.5 * phi * CLHEP::deg);
834 }
835 for (int i = 0; i < ndiv; ++i) {
836 const string physicalName = "physicalRib_" + to_string(id) + " " + to_string(i);
837 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
838 physicalName.c_str(), m_logicalCDC, false, id);
839 rot.rotateZ(phi * CLHEP::deg);
840 arm.rotateZ(phi * CLHEP::deg);
841 }
842
843 }
844
845 //
846 // Construct rib2s.
847 //
848 for (const auto& rib2 : geo.getRib2s()) {
849
850 const int id = rib2.getId();
851 const double length = rib2.getLength();
852 const double width = rib2.getWidth();
853 const double thick = rib2.getThick();
854 const double width2 = rib2.getWidth2();
855 const double thick2 = rib2.getThick2();
856 const double rotx = rib2.getRotX();
857 const double roty = rib2.getRotY();
858 const double rotz = rib2.getRotZ();
859 const double x = rib2.getX();
860 const double y = rib2.getY();
861 const double z = rib2.getZ();
862 const int ndiv = rib2.getNDiv();
863
864 const string solidName = "solidRib2" + to_string(id);
865 const string logicalName = "logicalRib2" + to_string(id);
866 G4Trd* trdShape = new G4Trd(solidName,
867 0.5 * thick * CLHEP::cm,
868 0.5 * thick2 * CLHEP::cm,
869 0.5 * width * CLHEP::cm,
870 0.5 * width2 * CLHEP::cm,
871 0.5 * length * CLHEP::cm);
872
873 G4LogicalVolume* logicalV = new G4LogicalVolume(trdShape, medAluminum, logicalName, 0, 0, 0);
874
875 if (id > 0)
876 logicalV = new G4LogicalVolume(trdShape, medCopper, logicalName, 0, 0, 0);
877
878 logicalV->SetVisAttributes(m_VisAttributes.back());
879
880 const double phi = 360.0 / ndiv;
881
882 G4RotationMatrix rot = G4RotationMatrix();
883 G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - thick * CLHEP::cm / 2.0);
884
885 rot.rotateX(rotx);
886 rot.rotateY(roty);
887 rot.rotateZ(rotz);
888 for (int i = 0; i < ndiv; ++i) {
889 const string physicalName = "physicalRib2_" + to_string(id) + " " + to_string(i);
890 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
891 physicalName.c_str(), m_logicalCDC, false, id);
892 rot.rotateZ(phi * CLHEP::deg);
893 arm.rotateZ(phi * CLHEP::deg);
894 }
895
896 }
897
898 //
899 // Construct rib3s.
900 //
901 for (const auto& rib3 : geo.getRib3s()) {
902
903 const int id = rib3.getId();
904 const double length = rib3.getLength();
905 const double width = rib3.getWidth();
906 const double thick = rib3.getThick();
907 const double r = rib3.getR();
908 const double x = rib3.getX();
909 const double y = rib3.getY();
910 const double z = rib3.getZ();
911 const double rx = rib3.getRx();
912 const double ry = rib3.getRy();
913 const double rz = rib3.getRz();
914 const int offset = rib3.getOffset();
915 const int ndiv = rib3.getNDiv();
916
917 const string logicalName = "logicalRib3" + to_string(id);
918 G4VSolid* boxShape = new G4Box("Block",
919 0.5 * length * CLHEP::cm,
920 0.5 * width * CLHEP::cm,
921 0.5 * thick * CLHEP::cm);
922 G4VSolid* tubeShape = new G4Tubs("Hole",
923 0.,
924 r * CLHEP::cm,
925 width * CLHEP::cm,
926 0. * CLHEP::deg,
927 360. * CLHEP::deg);
928
929 G4RotationMatrix rotsub = G4RotationMatrix();
930 rotsub.rotateX(90. * CLHEP::deg);
931 G4ThreeVector trnsub(rx * CLHEP::cm - x * CLHEP::cm, ry * CLHEP::cm - y * CLHEP::cm,
932 rz * CLHEP::cm - z * CLHEP::cm + 0.5 * thick * CLHEP::cm);
933 G4VSolid* coolingBlock = new G4SubtractionSolid("Block-Hole",
934 boxShape,
935 tubeShape,
936 G4Transform3D(rotsub,
937 trnsub));
938
939 G4LogicalVolume* logicalV = new G4LogicalVolume(coolingBlock, medCopper, logicalName, 0, 0, 0);
940
941 logicalV->SetVisAttributes(m_VisAttributes.back());
942
943 const double phi = 360.0 / ndiv;
944
945 G4RotationMatrix rot = G4RotationMatrix();
946 G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - thick * CLHEP::cm / 2.0);
947
948 if (offset) {
949 rot.rotateZ(0.5 * phi * CLHEP::deg);
950 arm.rotateZ(0.5 * phi * CLHEP::deg);
951 }
952 for (int i = 0; i < ndiv; ++i) {
953 const string physicalName = "physicalRib3_" + to_string(id) + " " + to_string(i);
954 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
955 physicalName.c_str(), m_logicalCDC, false, id);
956 rot.rotateZ(phi * CLHEP::deg);
957 arm.rotateZ(phi * CLHEP::deg);
958 }
959
960 }
961
962 //
963 // Construct rib4s.
964 //
965 for (const auto& rib4 : geo.getRib4s()) {
966
967 const int id = rib4.getId();
968 const double length = rib4.getLength();
969 const double width = rib4.getWidth();
970 const double thick = rib4.getThick();
971 const double length2 = rib4.getLength2();
972 const double width2 = rib4.getWidth2();
973 const double thick2 = rib4.getThick2();
974 const double x = rib4.getX();
975 const double y = rib4.getY();
976 const double z = rib4.getZ();
977 const double x2 = rib4.getX2();
978 const double y2 = rib4.getY2();
979 const double z2 = rib4.getZ2();
980 const int offset = rib4.getOffset();
981 const int ndiv = rib4.getNDiv();
982
983 const string logicalName = "logicalRib4" + to_string(id);
984 G4VSolid* baseShape = new G4Box("Base",
985 0.5 * length * CLHEP::cm,
986 0.5 * width * CLHEP::cm,
987 0.5 * thick * CLHEP::cm);
988 G4VSolid* sqShape = new G4Box("Sq",
989 0.5 * length2 * CLHEP::cm,
990 0.5 * width2 * CLHEP::cm,
991 0.5 * thick2 * CLHEP::cm);
992
993 G4RotationMatrix rotsub = G4RotationMatrix();
994 double dzc = (z2 - thick2 / 2.) - (z - thick / 2.);
995 G4ThreeVector trnsub(x2 * CLHEP::cm - x * CLHEP::cm,
996 y2 * CLHEP::cm - y * CLHEP::cm,
997 dzc * CLHEP::cm);
998 G4VSolid* sqHoleBase = new G4SubtractionSolid("Box-Sq",
999 baseShape,
1000 sqShape,
1001 G4Transform3D(rotsub,
1002 trnsub)
1003 );
1004
1005 G4LogicalVolume* logicalV = new G4LogicalVolume(sqHoleBase, medCopper, logicalName, 0, 0, 0);
1006 if (id < 19) {
1007 logicalV = new G4LogicalVolume(sqHoleBase, medNEMA_G10_Plate, logicalName, 0, 0, 0);
1008 BkgSensitiveDetector* sensitiveDetector =
1009 new BkgSensitiveDetector("CDC", 2000 + id);
1010 logicalV->SetSensitiveDetector(sensitiveDetector);
1011 m_BkgSensitiveRib4.push_back(sensitiveDetector);
1012 }
1013
1014 logicalV->SetVisAttributes(m_VisAttributes.back());
1015
1016 const double phi = 360.0 / ndiv;
1017
1018 G4RotationMatrix rot = G4RotationMatrix();
1019 G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - thick * CLHEP::cm / 2.0);
1020
1021 if (offset) {
1022 rot.rotateZ(0.5 * phi * CLHEP::deg);
1023 arm.rotateZ(0.5 * phi * CLHEP::deg);
1024 }
1025 for (int i = 0; i < ndiv; ++i) {
1026 const string physicalName = "physicalRib4_" + to_string(id) + " " + to_string(i);
1027 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1028 physicalName.c_str(), m_logicalCDC, false, id);
1029 rot.rotateZ(phi * CLHEP::deg);
1030 arm.rotateZ(phi * CLHEP::deg);
1031 }
1032
1033 }
1034 //
1035 // Construct rib5s.
1036 //
1037 for (const auto& rib5 : geo.getRib5s()) {
1038
1039 const int id = rib5.getId();
1040 const double dr = rib5.getDr();
1041 const double dz = rib5.getDz();
1042 const double width = rib5.getWidth();
1043 const double thick = rib5.getThick();
1044 const double rin = rib5.getRin();
1045 const double x = rib5.getX();
1046 const double y = rib5.getY();
1047 const double z = rib5.getZ();
1048 const double rotx = rib5.getRotx();
1049 const double roty = rib5.getRoty();
1050 const double rotz = rib5.getRotz();
1051 const int offset = rib5.getOffset();
1052 const int ndiv = rib5.getNDiv();
1053
1054 const string solidName = "solidRib5" + to_string(id);
1055 const string logicalName = "logicalRib5" + to_string(id);
1056
1057 const double rmax = rin + thick;
1058 const double rmin = rin;
1059 const double dphi = 2. * atan2(dz, dr);
1060 const double ddphi = thick * tan(dphi) / rin;
1061 const double ddphi2 = width / 2. * width / 2. / (x + dr) / rin;
1062 const double cphi = dphi - ddphi - ddphi2;
1063 G4Tubs* tubeShape = new G4Tubs(solidName,
1064 rmin * CLHEP::cm,
1065 rmax * CLHEP::cm,
1066 0.5 * width * CLHEP::cm,
1067 0.,
1068 cphi);
1069
1070 G4LogicalVolume* logicalV = new G4LogicalVolume(tubeShape, medAluminum, logicalName, 0, 0, 0);
1071
1072 logicalV->SetVisAttributes(m_VisAttributes.back());
1073
1074 const double phi = 360.0 / ndiv;
1075
1076 G4RotationMatrix rot = G4RotationMatrix();
1077
1078 //G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - thick * CLHEP::cm / 2.0);
1079 G4ThreeVector arm(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - rin * CLHEP::cm - thick * CLHEP::cm);
1080 rot.rotateX(rotx);
1081 rot.rotateY(roty);
1082 rot.rotateZ(rotz);
1083 if (offset) {
1084 rot.rotateZ(0.5 * phi * CLHEP::deg);
1085 arm.rotateZ(0.5 * phi * CLHEP::deg);
1086 }
1087 for (int i = 0; i < ndiv; ++i) {
1088 const string physicalName = "physicalRib5_" + to_string(id) + " " + to_string(i);
1089 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1090 physicalName.c_str(), m_logicalCDC, false, id);
1091 rot.rotateZ(phi * CLHEP::deg);
1092 arm.rotateZ(phi * CLHEP::deg);
1093 }
1094
1095 }
1096
1097 //Create B-field mapper (here tentatively)
1098 createMapper(topVolume);
1099 }
1100
1101
1103 {
1104
1105 G4Material* C2H4 = geometry::Materials::get("G4_POLYETHYLENE");
1106 G4Material* elB = geometry::Materials::get("G4_B");
1107
1108 // 5% borated polyethylene = SWX201
1109 // http://www.deqtech.com/Shieldwerx/Products/swx201hd.htm
1110 G4Material* boratedpoly05 = new G4Material("BoratedPoly05", 1.06 * CLHEP::g / CLHEP::cm3, 2);
1111 boratedpoly05->AddMaterial(elB, 0.05);
1112 boratedpoly05->AddMaterial(C2H4, 0.95);
1113 // 30% borated polyethylene = SWX210
1114 G4Material* boratedpoly30 = new G4Material("BoratedPoly30", 1.19 * CLHEP::g / CLHEP::cm3, 2);
1115 boratedpoly30->AddMaterial(elB, 0.30);
1116 boratedpoly30->AddMaterial(C2H4, 0.70);
1117
1118 G4Material* shieldMat = C2H4;
1119
1120 const int nShields = content.getNumberNodes("Shields/Shield");
1121
1122 for (int iShield = 0; iShield < nShields; ++iShield) {
1123 GearDir shieldContent(content);
1124 shieldContent.append((format("/Shields/Shield[%1%]/") % (iShield + 1)).str());
1125 const string sShieldID = shieldContent.getString("@id");
1126 const int shieldID = atoi(sShieldID.c_str());
1127 // const string shieldName = shieldContent.getString("Name");
1128 const double shieldInnerR1 = shieldContent.getLength("InnerR1");
1129 const double shieldInnerR2 = shieldContent.getLength("InnerR2");
1130 const double shieldOuterR1 = shieldContent.getLength("OuterR1");
1131 const double shieldOuterR2 = shieldContent.getLength("OuterR2");
1132 const double shieldThick = shieldContent.getLength("Thickness");
1133 const double shieldPosZ = shieldContent.getLength("PosZ");
1134
1135 G4Cons* shieldConsShape = new G4Cons((format("solidShield%1%") % shieldID).str().c_str(),
1136 shieldInnerR1 * CLHEP::cm, shieldOuterR1 * CLHEP::cm,
1137 shieldInnerR2 * CLHEP::cm, shieldOuterR2 * CLHEP::cm,
1138 shieldThick * CLHEP::cm / 2.0,
1139 0.*CLHEP::deg, 360.*CLHEP::deg);
1140
1141 G4LogicalVolume* shieldCons = new G4LogicalVolume(shieldConsShape, shieldMat, (format("logicalShield%1%") % shieldID).str().c_str(),
1142 0, 0, 0);
1143 shieldCons->SetVisAttributes(m_VisAttributes.back());
1144 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (shieldPosZ - shieldThick / 2.0) * CLHEP::cm), shieldCons,
1145 (format("physicalShield%1%") % shieldID).str().c_str(), m_logicalCDC, false, 0);
1146
1147 }
1148
1149 }
1150
1151
1153 {
1154
1155 G4Material* C2H4 = geometry::Materials::get("G4_POLYETHYLENE");
1156 G4Material* shieldMat = C2H4;
1157
1158 for (const auto& shield : geom.getNeutronShields()) {
1159 const int shieldID = shield.getId();
1160 const double shieldInnerR1 = shield.getRmin1();
1161 const double shieldInnerR2 = shield.getRmin2();
1162 const double shieldOuterR1 = shield.getRmax1();
1163 const double shieldOuterR2 = shield.getRmax2();
1164 const double shieldThick = shield.getThick();
1165 const double shieldPosZ = shield.getZ();
1166
1167 G4Cons* shieldConsShape = new G4Cons("solidShield" + to_string(shieldID),
1168 shieldInnerR1 * CLHEP::cm, shieldOuterR1 * CLHEP::cm,
1169 shieldInnerR2 * CLHEP::cm, shieldOuterR2 * CLHEP::cm,
1170 shieldThick * CLHEP::cm / 2.0,
1171 0.*CLHEP::deg, 360.*CLHEP::deg);
1172
1173 G4LogicalVolume* shieldCons = new G4LogicalVolume(shieldConsShape, shieldMat, "logicalShield" + to_string(shieldID),
1174 0, 0, 0);
1175 shieldCons->SetVisAttributes(m_VisAttributes.back());
1176 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, (shieldPosZ - shieldThick / 2.0) * CLHEP::cm), shieldCons,
1177 "physicalShield" + to_string(shieldID), m_logicalCDC, false, 0);
1178
1179 }
1180
1181 }
1182
1184 {
1185 string Aluminum = content.getString("Aluminum");
1186 G4Material* medAluminum = geometry::Materials::get(Aluminum);
1187 G4Material* medNEMA_G10_Plate = geometry::Materials::get("NEMA_G10_Plate");
1188 G4double density = 1.000 * CLHEP::g / CLHEP::cm3;
1189 G4double a = 1.01 * CLHEP::g / CLHEP::mole;
1190 G4Element* elH = new G4Element("Hydrogen", "H", 1., a);
1191 a = 16.00 * CLHEP::g / CLHEP::mole;
1192 G4Element* elO = new G4Element("Oxygen", "O", 8., a);
1193 G4Material* medH2O = new G4Material("Water", density, 2);
1194 medH2O->AddElement(elH, 2);
1195 medH2O->AddElement(elO, 1);
1196 G4Material* medCopper = geometry::Materials::get("Cu");
1197 G4Material* medLV = geometry::Materials::get("CDCLVCable");
1198 G4Material* medFiber = geometry::Materials::get("CDCOpticalFiber");
1199 G4Material* medCAT7 = geometry::Materials::get("CDCCAT7");
1200 G4Material* medTRG = geometry::Materials::get("CDCOpticalFiberTRG");
1201 G4Material* medHV = geometry::Materials::get("CDCHVCable");
1202
1203 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(0., 1., 0.)));
1204 const int nCover = content.getNumberNodes("Covers/Cover");
1205 for (int iCover = 0; iCover < nCover; ++iCover) {
1206 GearDir coverContent(content);
1207 coverContent.append((format("/Covers/Cover[%1%]/") % (iCover + 1)).str());
1208 const string scoverID = coverContent.getString("@id");
1209 const int coverID = atoi(scoverID.c_str());
1210 const string coverName = coverContent.getString("Name");
1211 const double coverInnerR1 = coverContent.getLength("InnerR1");
1212 const double coverInnerR2 = coverContent.getLength("InnerR2");
1213 const double coverOuterR1 = coverContent.getLength("OuterR1");
1214 const double coverOuterR2 = coverContent.getLength("OuterR2");
1215 const double coverThick = coverContent.getLength("Thickness");
1216 const double coverPosZ = coverContent.getLength("PosZ");
1217
1218 const double rmin1 = coverInnerR1;
1219 const double rmax1 = coverOuterR1;
1220 const double rmin2 = coverInnerR2;
1221 const double rmax2 = coverOuterR2;
1222
1223 /*
1224 if (coverID == 7 || coverID == 10) {
1225 createCone(rmin1, rmax1, rmin2, rmax2, coverThick, coverPosZ, coverID, medAluminum, coverName);
1226 } else {
1227 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medAluminum, coverName);
1228
1229 }*/
1230 // ID dependent material definition
1231 if (coverID < 23) {
1232 if (coverID == 7 || coverID == 10) {// cones
1233 createCone(rmin1, rmax1, rmin2, rmax2, coverThick, coverPosZ, coverID, medAluminum, coverName);
1234 } else {// covers
1235 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medAluminum, coverName);
1236 }
1237 }
1238 if (coverID > 22 && coverID < 29)// cooling plate
1239 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medCopper, coverName);
1240 if (coverID > 28 && coverID < 35)// cooling Pipe
1241 createTorus(rmin1, rmax1, coverThick, coverPosZ, coverID, medCopper, coverName);
1242 if (coverID > 34 && coverID < 41)// cooling water
1243 createTorus(rmin1, rmax1, coverThick, coverPosZ, coverID, medH2O, coverName);
1244 if (coverID == 45 || coverID == 46)
1245 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medLV, coverName);
1246 if (coverID == 47 || coverID == 48)
1247 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medFiber, coverName);
1248 if (coverID == 49 || coverID == 50)
1249 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medCAT7, coverName);
1250 if (coverID == 51 || coverID == 52)
1251 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medTRG, coverName);
1252 if (coverID == 53)
1253 createTube(rmin1, rmax1, coverThick, coverPosZ, coverID, medHV, coverName);
1254 }
1255
1256 const int nCover2 = content.getNumberNodes("Covers/Cover2");
1257 for (int iCover2 = 0; iCover2 < nCover2; ++iCover2) {
1258 GearDir cover2Content(content);
1259 cover2Content.append((format("/Cover2s/Cover2[%1%]/") % (iCover2 + 1)).str());
1260 const string scover2ID = cover2Content.getString("@id");
1261 const int cover2ID = atoi(scover2ID.c_str());
1262 const string cover2Name = cover2Content.getString("Name");
1263 const double cover2InnerR = cover2Content.getLength("InnerR");
1264 const double cover2OuterR = cover2Content.getLength("OuterR");
1265 const double cover2StartPhi = cover2Content.getLength("StartPhi");
1266 const double cover2DeltaPhi = cover2Content.getLength("DeltaPhi");
1267 const double cover2Thick = cover2Content.getLength("Thickness");
1268 const double cover2PosZ = cover2Content.getLength("PosZ");
1269
1270 if (cover2ID < 11)
1271 createTube2(cover2InnerR, cover2OuterR, cover2StartPhi, cover2DeltaPhi, cover2Thick, cover2PosZ, cover2ID, medHV, cover2Name);
1272 if (cover2ID > 10 && cover2ID < 14)
1273 createTube2(cover2InnerR, cover2OuterR, cover2StartPhi, cover2DeltaPhi, cover2Thick, cover2PosZ, cover2ID, medFiber, cover2Name);
1274 if (cover2ID > 13 && cover2ID < 23)
1275 createTube2(cover2InnerR, cover2OuterR, cover2StartPhi, cover2DeltaPhi, cover2Thick, cover2PosZ, cover2ID, medCAT7, cover2Name);
1276 if (cover2ID > 22 && cover2ID < 29)
1277 createTube2(cover2InnerR, cover2OuterR, cover2StartPhi, cover2DeltaPhi, cover2Thick, cover2PosZ, cover2ID, medTRG, cover2Name);
1278 }
1279
1280 const int nRibs = content.getNumberNodes("Covers/Rib");
1281 for (int iRib = 0; iRib < nRibs; ++iRib) {
1282 GearDir ribContent(content);
1283 ribContent.append((format("/Covers/Rib[%1%]/") % (iRib + 1)).str());
1284 const string sribID = ribContent.getString("@id");
1285 const int ribID = atoi(sribID.c_str());
1286 // const string ribName = ribContent.getString("Name");
1287 const double length = ribContent.getLength("Length");
1288 const double width = ribContent.getLength("Width");
1289 const double thick = ribContent.getLength("Thickness");
1290 const double rotX = ribContent.getLength("RotX");
1291 const double rotY = ribContent.getLength("RotY");
1292 const double rotZ = ribContent.getLength("RotZ");
1293 const double cX = ribContent.getLength("PosX");
1294 const double cY = ribContent.getLength("PosY");
1295 const double cZ = ribContent.getLength("PosZ");
1296 const int offset = atoi((ribContent.getString("Offset")).c_str());
1297 const int number = atoi((ribContent.getString("NDiv")).c_str());
1298
1299 const string solidName = "solidRib" + to_string(ribID);
1300 const string logicalName = "logicalRib" + to_string(ribID);
1301 G4Box* boxShape = new G4Box(solidName, 0.5 * length * CLHEP::cm,
1302 0.5 * width * CLHEP::cm,
1303 0.5 * thick * CLHEP::cm);
1304 const double rmax = 0.5 * length;
1305 const double rmin = max((rmax - thick), 0.);
1306 G4Tubs* tubeShape = new G4Tubs(solidName,
1307 rmin * CLHEP::cm,
1308 rmax * CLHEP::cm,
1309 0.5 * width * CLHEP::cm,
1310 0.,
1311 360. * CLHEP::deg);
1312
1313 //G4LogicalVolume* logicalV = new G4LogicalVolume(boxShape, medAluminum,
1314 // logicalName, 0, 0, 0);
1315 // ID dependent material definition Aluminum is default
1316 G4LogicalVolume* logicalV = new G4LogicalVolume(boxShape, medAluminum, logicalName, 0, 0, 0);
1317 if (ribID > 39 && ribID < 78) // Cu box
1318 logicalV = new G4LogicalVolume(boxShape, medCopper, logicalName, 0, 0, 0);
1319 if ((ribID > 77 && ribID < 94) || (ribID > 131 && ribID < 146)) // G10 box
1320 logicalV = new G4LogicalVolume(boxShape, medNEMA_G10_Plate, logicalName, 0, 0, 0);
1321 if (ribID > 93 && ribID < 110) // Cu tube
1322 logicalV = new G4LogicalVolume(tubeShape, medCopper, logicalName, 0, 0, 0);
1323 if (ribID > 109 && ribID < 126) // H2O tube (rmin = 0)
1324 logicalV = new G4LogicalVolume(tubeShape, medH2O, logicalName, 0, 0, 0);
1325 if (ribID > 127 && ribID < 132) // HV bundle
1326 logicalV = new G4LogicalVolume(boxShape, medHV, logicalName, 0, 0, 0);
1327 /*if( ribID > 145 && ribID < 149 )// Fiber box
1328 logicalV = new G4LogicalVolume(boxShape, medFiber, logicalName, 0, 0, 0);
1329 if( ribID > 148 && ribID < 158 )// Fiber box
1330 logicalV = new G4LogicalVolume(boxShape, medCAT7, logicalName, 0, 0, 0);
1331 if( ribID > 157 && ribID < 164 )// Fiber box
1332 logicalV = new G4LogicalVolume(boxShape, medTRG, logicalName, 0, 0, 0);*/
1333
1334 logicalV->SetVisAttributes(m_VisAttributes.back());
1335
1336 const double phi = 360.0 / number;
1337
1338 G4RotationMatrix rot = G4RotationMatrix();
1339
1340 double dz = thick;
1341 if (ribID > 93 && ribID < 126) dz = 0;
1342 G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - dz * CLHEP::cm / 2.0);
1343
1344 rot.rotateX(rotX);
1345 rot.rotateY(rotY);
1346 rot.rotateZ(rotZ);
1347 if (offset) {
1348 rot.rotateZ(0.5 * phi * CLHEP::deg);
1349 arm.rotateZ(0.5 * phi * CLHEP::deg);
1350 }
1351 for (int i = 0; i < number; ++i) {
1352 const string physicalName = "physicalRib_" + to_string(ribID) + " " + to_string(i);
1353 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1354 physicalName.c_str(), m_logicalCDC, false, ribID);
1355 rot.rotateZ(phi * CLHEP::deg);
1356 arm.rotateZ(phi * CLHEP::deg);
1357 }
1358
1359 }// rib
1360
1361 const int nRib2s = content.getNumberNodes("Covers/Rib2");
1362 for (int iRib2 = 0; iRib2 < nRib2s; ++iRib2) {
1363 GearDir rib2Content(content);
1364 rib2Content.append((format("/Covers/Rib2[%1%]/") % (iRib2 + 1)).str());
1365 const string srib2ID = rib2Content.getString("@id");
1366 const int rib2ID = atoi(srib2ID.c_str());
1367 // const string rib2Name = rib2Content.getString("Name");
1368 const double length = rib2Content.getLength("Length");
1369 const double width = rib2Content.getLength("Width");
1370 const double thick = rib2Content.getLength("Thickness");
1371 const double width2 = rib2Content.getLength("Width2");
1372 const double thick2 = rib2Content.getLength("Thickness2");
1373 const double rotX = rib2Content.getLength("RotX");
1374 const double rotY = rib2Content.getLength("RotY");
1375 const double rotZ = rib2Content.getLength("RotZ");
1376 const double cX = rib2Content.getLength("PosX");
1377 const double cY = rib2Content.getLength("PosY");
1378 const double cZ = rib2Content.getLength("PosZ");
1379 const int number = atoi((rib2Content.getString("NDiv")).c_str());
1380
1381 const string solidName = "solidRib2" + to_string(rib2ID);
1382 const string logicalName = "logicalRib2" + to_string(rib2ID);
1383 G4Trd* trdShape = new G4Trd(solidName,
1384 0.5 * thick * CLHEP::cm,
1385 0.5 * thick2 * CLHEP::cm,
1386 0.5 * width * CLHEP::cm,
1387 0.5 * width2 * CLHEP::cm,
1388 0.5 * length * CLHEP::cm);
1389
1390 G4LogicalVolume* logicalV = new G4LogicalVolume(trdShape, medAluminum, logicalName, 0, 0, 0);
1391 if (rib2ID > 0)
1392 logicalV = new G4LogicalVolume(trdShape, medCopper, logicalName, 0, 0, 0);
1393
1394 logicalV->SetVisAttributes(m_VisAttributes.back());
1395
1396 const double phi = 360.0 / number;
1397
1398 G4RotationMatrix rot = G4RotationMatrix();
1399 G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - thick * CLHEP::cm / 2.0);
1400
1401 rot.rotateX(rotX);
1402 rot.rotateY(rotY);
1403 rot.rotateZ(rotZ);
1404 for (int i = 0; i < number; ++i) {
1405 const string physicalName = "physicalRib2_" + to_string(rib2ID) + " " + to_string(i);
1406 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1407 physicalName.c_str(), m_logicalCDC, false, rib2ID);
1408 rot.rotateZ(phi * CLHEP::deg);
1409 arm.rotateZ(phi * CLHEP::deg);
1410 }
1411
1412 }// rib2
1413
1414 const int nRib3s = content.getNumberNodes("Covers/Rib3");
1415 for (int iRib3 = 0; iRib3 < nRib3s; ++iRib3) {
1416 GearDir rib3Content(content);
1417 rib3Content.append((format("/Covers/Rib3[%1%]/") % (iRib3 + 1)).str());
1418 const string srib3ID = rib3Content.getString("@id");
1419 const int rib3ID = atoi(srib3ID.c_str());
1420 // const string rib3Name = rib3Content.getString("Name");
1421 const double length = rib3Content.getLength("Length");
1422 const double width = rib3Content.getLength("Width");
1423 const double thick = rib3Content.getLength("Thickness");
1424 const double r = rib3Content.getLength("HoleR");
1425 const double cX = rib3Content.getLength("PosX");
1426 const double cY = rib3Content.getLength("PosY");
1427 const double cZ = rib3Content.getLength("PosZ");
1428 const double hX = rib3Content.getLength("HoleX");
1429 const double hY = rib3Content.getLength("HoleY");
1430 const double hZ = rib3Content.getLength("HoleZ");
1431 const int offset = atoi((rib3Content.getString("Offset")).c_str());
1432 const int number = atoi((rib3Content.getString("NDiv")).c_str());
1433
1434 const string logicalName = "logicalRib3" + to_string(rib3ID);
1435 G4VSolid* boxShape = new G4Box("Block",
1436 0.5 * length * CLHEP::cm,
1437 0.5 * width * CLHEP::cm,
1438 0.5 * thick * CLHEP::cm);
1439 G4VSolid* tubeShape = new G4Tubs("Hole",
1440 0.,
1441 r * CLHEP::cm,
1442 length * CLHEP::cm,
1443 0.,
1444 360. * CLHEP::deg);
1445 G4RotationMatrix rotsub = G4RotationMatrix();
1446 G4ThreeVector trnsub(cX * CLHEP::cm - hX * CLHEP::cm, cY * CLHEP::cm - hY * CLHEP::cm,
1447 cZ * CLHEP::cm - hZ * CLHEP::cm + 0.5 * thick * CLHEP::cm);
1448 G4VSolid* coolingBlock = new G4SubtractionSolid("Block-Hole",
1449 boxShape,
1450 tubeShape,
1451 G4Transform3D(rotsub,
1452 trnsub));
1453
1454 G4LogicalVolume* logicalV = new G4LogicalVolume(coolingBlock, medCopper, logicalName, 0, 0, 0);
1455
1456 logicalV->SetVisAttributes(m_VisAttributes.back());
1457
1458 const double phi = 360.0 / number;
1459
1460 G4RotationMatrix rot = G4RotationMatrix();
1461 G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - thick * CLHEP::cm / 2.0);
1462
1463 if (offset) {
1464 rot.rotateZ(0.5 * phi * CLHEP::deg);
1465 arm.rotateZ(0.5 * phi * CLHEP::deg);
1466 }
1467 for (int i = 0; i < number; ++i) {
1468 const string physicalName = "physicalRib3_" + to_string(rib3ID) + " " + to_string(i);
1469 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1470 physicalName.c_str(), m_logicalCDC, false, rib3ID);
1471 rot.rotateZ(phi * CLHEP::deg);
1472 arm.rotateZ(phi * CLHEP::deg);
1473 }
1474
1475 }// rib3
1476
1477 const int nRib4s = content.getNumberNodes("Covers/Rib4");
1478 for (int iRib4 = 0; iRib4 < nRib4s; ++iRib4) {
1479 GearDir rib4Content(content);
1480 rib4Content.append((format("/Covers/Rib4[%1%]/") % (iRib4 + 1)).str());
1481 const string srib4ID = rib4Content.getString("@id");
1482 const int rib4ID = atoi(srib4ID.c_str());
1483 // const string rib4Name = rib4Content.getString("Name");
1484 const double length = rib4Content.getLength("Length");
1485 const double width = rib4Content.getLength("Width");
1486 const double thick = rib4Content.getLength("Thickness");
1487 const double length2 = rib4Content.getLength("Length2");
1488 const double width2 = rib4Content.getLength("Width2");
1489 const double thick2 = rib4Content.getLength("Thickness2");
1490 const double cX = rib4Content.getLength("PosX");
1491 const double cY = rib4Content.getLength("PosY");
1492 const double cZ = rib4Content.getLength("PosZ");
1493 const double hX = rib4Content.getLength("HoleX");
1494 const double hY = rib4Content.getLength("HoleY");
1495 const double hZ = rib4Content.getLength("HoleZ");
1496 const int offset = atoi((rib4Content.getString("Offset")).c_str());
1497 const int number = atoi((rib4Content.getString("NDiv")).c_str());
1498
1499 const string logicalName = "logicalRib4" + to_string(rib4ID);
1500 G4VSolid* baseShape = new G4Box("Base",
1501 0.5 * length * CLHEP::cm,
1502 0.5 * width * CLHEP::cm,
1503 0.5 * thick * CLHEP::cm);
1504 G4VSolid* sqShape = new G4Box("Sq",
1505 0.5 * length2 * CLHEP::cm,
1506 0.5 * width2 * CLHEP::cm,
1507 0.5 * thick2 * CLHEP::cm);
1508 G4RotationMatrix rotsub = G4RotationMatrix();
1509 double dzc = (hZ - thick2 / 2.) - (cZ - thick / 2.);
1510 G4ThreeVector trnsub(hX * CLHEP::cm - cX * CLHEP::cm,
1511 hY * CLHEP::cm - cY * CLHEP::cm,
1512 dzc * CLHEP::cm);
1513 G4VSolid* sqHoleBase = new G4SubtractionSolid("Base-Sq",
1514 baseShape,
1515 sqShape,
1516 G4Transform3D(rotsub,
1517 trnsub)
1518 );
1519
1520 G4LogicalVolume* logicalV = new G4LogicalVolume(sqHoleBase, medCopper, logicalName, 0, 0, 0);
1521 if (rib4ID < 19)
1522 logicalV = new G4LogicalVolume(sqHoleBase, medNEMA_G10_Plate, logicalName, 0, 0, 0);
1523
1524 logicalV->SetVisAttributes(m_VisAttributes.back());
1525
1526 const double phi = 360.0 / number;
1527
1528 G4RotationMatrix rot = G4RotationMatrix();
1529 G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - thick * CLHEP::cm / 2.0);
1530
1531 if (offset) {
1532 rot.rotateZ(0.5 * phi * CLHEP::deg);
1533 arm.rotateZ(0.5 * phi * CLHEP::deg);
1534 }
1535 for (int i = 0; i < number; ++i) {
1536 const string physicalName = "physicalRib4_" + to_string(rib4ID) + " " + to_string(i);
1537 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1538 physicalName.c_str(), m_logicalCDC, false, rib4ID);
1539 rot.rotateZ(phi * CLHEP::deg);
1540 arm.rotateZ(phi * CLHEP::deg);
1541 }
1542
1543 }// rib4
1544
1545 const int nRib5s = content.getNumberNodes("Covers/Rib5");
1546 for (int iRib5 = 0; iRib5 < nRib5s; ++iRib5) {
1547 GearDir rib5Content(content);
1548 rib5Content.append((format("/Covers/Rib5[%1%]/") % (iRib5 + 1)).str());
1549 const string srib5ID = rib5Content.getString("@id");
1550 const int rib5ID = atoi(srib5ID.c_str());
1551 // const string rib5Name = rib5Content.getString("Name");
1552 const double dr = rib5Content.getLength("DeltaR");
1553 const double dz = rib5Content.getLength("DeltaZ");
1554 const double width = rib5Content.getLength("Width");
1555 const double thick = rib5Content.getLength("Thickness");
1556 const double rin = rib5Content.getLength("Rin");
1557 const double rotX = rib5Content.getLength("RotX");
1558 const double rotY = rib5Content.getLength("RotY");
1559 const double rotZ = rib5Content.getLength("RotZ");
1560 const double cX = rib5Content.getLength("PosX");
1561 const double cY = rib5Content.getLength("PosY");
1562 const double cZ = rib5Content.getLength("PosZ");
1563 const int offset = atoi((rib5Content.getString("Offset")).c_str());
1564 const int number = atoi((rib5Content.getString("NDiv")).c_str());
1565
1566 const string solidName = "solidRib5" + to_string(rib5ID);
1567 const string logicalName = "logicalRib5" + to_string(rib5ID);
1568 const double rmax = rin + thick;
1569 const double rmin = rin;
1570 const double dphi = 2. * atan2(dz, dr);
1571 const double ddphi = thick * tan(dphi) / rin;
1572 const double ddphi2 = width / 2. * width / 2. / (cX + dr) / rin;
1573 const double cphi = dphi - ddphi - ddphi2;
1574 G4Tubs* tubeShape = new G4Tubs(solidName,
1575 rmin * CLHEP::cm,
1576 rmax * CLHEP::cm,
1577 0.5 * width * CLHEP::cm,
1578 0.,
1579 cphi);
1580
1581 G4LogicalVolume* logicalV = new G4LogicalVolume(tubeShape, medAluminum, logicalName, 0, 0, 0);
1582
1583 logicalV->SetVisAttributes(m_VisAttributes.back());
1584
1585 const double phi = 360.0 / number;
1586
1587 G4RotationMatrix rot = G4RotationMatrix();
1588
1589 //G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - thick * CLHEP::cm / 2.0);
1590 G4ThreeVector arm(cX * CLHEP::cm, cY * CLHEP::cm, cZ * CLHEP::cm - rin * CLHEP::cm - thick * CLHEP::cm);
1591
1592 rot.rotateX(rotX);
1593 rot.rotateY(rotY);
1594 rot.rotateZ(rotZ);
1595 if (offset) {
1596 rot.rotateZ(0.5 * phi * CLHEP::deg);
1597 arm.rotateZ(0.5 * phi * CLHEP::deg);
1598 }
1599 for (int i = 0; i < number; ++i) {
1600 const string physicalName = "physicalRib5_" + to_string(rib5ID) + " " + to_string(i);
1601 new G4PVPlacement(G4Transform3D(rot, arm), logicalV,
1602 physicalName.c_str(), m_logicalCDC, false, rib5ID);
1603 rot.rotateZ(phi * CLHEP::deg);
1604 arm.rotateZ(phi * CLHEP::deg);
1605 }
1606 }//rib5
1607
1608 }
1609
1610
1612 {
1613 G4Material* medAl = geometry::Materials::get("Al");
1614 G4double density = 1.000 * CLHEP::g / CLHEP::cm3;
1615 G4double a = 1.01 * CLHEP::g / CLHEP::mole;
1616 G4Element* elH = new G4Element("Hydrogen", "H", 1., a);
1617 a = 16.00 * CLHEP::g / CLHEP::mole;
1618 G4Element* elO = new G4Element("Oxygen", "O", 8., a);
1619 G4Material* medH2O = new G4Material("water", density, 2);
1620 medH2O->AddElement(elH, 2);
1621 medH2O->AddElement(elO, 1);
1622 G4Material* medCu = geometry::Materials::get("Cu");
1623 G4Material* medLV = geometry::Materials::get("CDCLVCable");
1624 G4Material* medFiber = geometry::Materials::get("CDCOpticalFiber");
1625 G4Material* medCAT7 = geometry::Materials::get("CDCCAT7");
1626 G4Material* medTRG = geometry::Materials::get("CDCOpticalFiberTRG");
1627 G4Material* medHV = geometry::Materials::get("CDCHVCable");
1628
1629 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(0., 1., 0.)));
1630 for (const auto& cover : geom.getCovers()) {
1631 const int coverID = cover.getId();
1632 const string coverName = "cover" + to_string(coverID);
1633 const double rmin1 = cover.getRmin1();
1634 const double rmin2 = cover.getRmin2();
1635 const double rmax1 = cover.getRmax1();
1636 const double rmax2 = cover.getRmax2();
1637 const double thick = cover.getThick();
1638 const double posZ = cover.getZ();
1639
1640 /*if (coverID == 7 || coverID == 10) {
1641 createCone(rmin1, rmax1, rmin2, rmax2, thick, posZ, coverID, medAl, coverName);
1642 } else {
1643 createTube(rmin1, rmax1, thick, posZ, coverID, medAl, coverName);
1644 }*/
1645 // ID dependent material definition
1646 if (coverID < 23) {
1647 if (coverID == 7 || coverID == 10) {
1648 createCone(rmin1, rmax1, rmin2, rmax2, thick, posZ, coverID, medAl, coverName);
1649 } else {
1650 createTube(rmin1, rmax1, thick, posZ, coverID, medAl, coverName);
1651 }
1652 }
1653 if (coverID > 22 && coverID < 29)
1654 createTube(rmin1, rmax1, thick, posZ, coverID, medCu, coverName);
1655 if (coverID > 28 && coverID < 35)
1656 createTorus(rmin1, rmax1, thick, posZ, coverID, medCu, coverName);
1657 if (coverID > 34 && coverID < 41)
1658 createTorus(rmin1, rmax1, thick, posZ, coverID, medH2O, coverName);
1659 if (coverID == 45 || coverID == 46)
1660 createTube(rmin1, rmax1, thick, posZ, coverID, medLV, coverName);
1661 if (coverID == 47 || coverID == 48)
1662 createTube(rmin1, rmax1, thick, posZ, coverID, medFiber, coverName);
1663 if (coverID == 49 || coverID == 50)
1664 createTube(rmin1, rmax1, thick, posZ, coverID, medCAT7, coverName);
1665 if (coverID == 51 || coverID == 52)
1666 createTube(rmin1, rmax1, thick, posZ, coverID, medTRG, coverName);
1667 if (coverID == 53)
1668 createTube(rmin1, rmax1, thick, posZ, coverID, medHV, coverName);
1669 }
1670 }
1671
1673 {
1674 G4Material* medHV = geometry::Materials::get("CDCHVCable");
1675 G4Material* medFiber = geometry::Materials::get("CDCOpticalFiber");
1676 G4Material* medCAT7 = geometry::Materials::get("CDCCAT7");
1677 G4Material* medTRG = geometry::Materials::get("CDCOpticalFiberTRG");
1678
1679 m_VisAttributes.push_back(new G4VisAttributes(true, G4Colour(0., 1., 0.)));
1680 for (const auto& cover2 : geom.getCover2s()) {
1681 const int cover2ID = cover2.getId();
1682 const string cover2Name = "cover2" + to_string(cover2ID);
1683 const double rmin = cover2.getRmin();
1684 const double rmax = cover2.getRmax();
1685 const double phis = cover2.getPhis();
1686 const double dphi = cover2.getDphi();
1687 const double thick = cover2.getThick();
1688 const double posZ = cover2.getZ();
1689
1690 if (cover2ID < 11)
1691 createTube2(rmin, rmax, phis, dphi, thick, posZ, cover2ID, medHV, cover2Name);
1692 if (cover2ID > 10 && cover2ID < 14)
1693 createTube2(rmin, rmax, phis, dphi, thick, posZ, cover2ID, medFiber, cover2Name);
1694 if (cover2ID > 13 && cover2ID < 23)
1695 createTube2(rmin, rmax, phis, dphi, thick, posZ, cover2ID, medCAT7, cover2Name);
1696 if (cover2ID > 22 && cover2ID < 29)
1697 createTube2(rmin, rmax, phis, dphi, thick, posZ, cover2ID, medTRG, cover2Name);
1698 }
1699 }
1700
1701 void GeoCDCCreator::createCone(const double rmin1, const double rmax1,
1702 const double rmin2, const double rmax2,
1703 const double thick, const double posZ,
1704 const int id, G4Material* med,
1705 const string& name)
1706 {
1707 const string solidName = "solid" + name;
1708 const string logicalName = "logical" + name;
1709 const string physicalName = "physical" + name;
1710 G4Cons* coverConeShape = new G4Cons(solidName.c_str(), rmin1 * CLHEP::cm, rmax1 * CLHEP::cm,
1711 rmin2 * CLHEP::cm, rmax2 * CLHEP::cm, thick * CLHEP::cm / 2.0, 0.*CLHEP::deg, 360.*CLHEP::deg);
1712 G4LogicalVolume* coverCone = new G4LogicalVolume(coverConeShape, med,
1713 logicalName.c_str(), 0, 0, 0);
1714 coverCone->SetVisAttributes(m_VisAttributes.back());
1715 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, posZ * CLHEP::cm - thick * CLHEP::cm / 2.0), coverCone,
1716 physicalName.c_str(), m_logicalCDC, false, id);
1717
1718 }
1719
1720 void GeoCDCCreator::createTube(const double rmin, const double rmax,
1721 const double thick, const double posZ,
1722 const int id, G4Material* med,
1723 const string& name)
1724 {
1725 const string solidName = "solid" + name;
1726 const string logicalName = "logical" + name;
1727 const string physicalName = "physical" + name;
1728 G4Tubs* solidV = new G4Tubs(solidName.c_str(),
1729 rmin * CLHEP::cm,
1730 rmax * CLHEP::cm,
1731 thick * CLHEP::cm / 2.0,
1732 0.*CLHEP::deg,
1733 360.*CLHEP::deg);
1734 G4LogicalVolume* logicalV = new G4LogicalVolume(solidV, med,
1735 logicalName.c_str(), 0, 0, 0);
1736 logicalV->SetVisAttributes(m_VisAttributes.back());
1737 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, posZ * CLHEP::cm - thick * CLHEP::cm / 2.0), logicalV,
1738 physicalName.c_str(), m_logicalCDC, false, id);
1739
1740 }
1741
1742 void GeoCDCCreator::createBox(const double length, const double height,
1743 const double thick, const double x,
1744 const double y, const double z,
1745 const int id, G4Material* med,
1746 const string& name)
1747 {
1748 const string solidName = (format("solid%1%%2%") % name % id).str();
1749 const string logicalName = (format("logical%1%%2%") % name % id).str();
1750 const string physicalName = (format("physical%1%%2%") % name % id).str();
1751 G4Box* boxShape = new G4Box(solidName.c_str(), 0.5 * length * CLHEP::cm,
1752 0.5 * height * CLHEP::cm,
1753 0.5 * thick * CLHEP::cm);
1754 G4LogicalVolume* logicalV = new G4LogicalVolume(boxShape, med,
1755 logicalName.c_str(), 0, 0, 0);
1756 logicalV->SetVisAttributes(m_VisAttributes.back());
1757 new G4PVPlacement(0, G4ThreeVector(x * CLHEP::cm, y * CLHEP::cm, z * CLHEP::cm - thick * CLHEP::cm / 2.0), logicalV,
1758 physicalName.c_str(), m_logicalCDC, false, id);
1759
1760 }
1761
1762 void GeoCDCCreator::createTorus(const double rmin1, const double rmax1,
1763 const double thick, const double posZ,
1764 const int id, G4Material* med,
1765 const string& name)
1766 {
1767 const string solidName = "solid" + name;
1768 const string logicalName = "logical" + name;
1769 const string physicalName = "physical" + name;
1770 const double rtor = (rmax1 + rmin1) / 2.;
1771 const double rmax = rmax1 - rtor;
1772 const double rmin = max((rmax - thick), 0.);
1773
1774 G4Torus* solidV = new G4Torus(solidName.c_str(),
1775 rmin * CLHEP::cm,
1776 rmax * CLHEP::cm,
1777 rtor * CLHEP::cm,
1778 0.*CLHEP::deg,
1779 360.*CLHEP::deg);
1780 G4LogicalVolume* logicalV = new G4LogicalVolume(solidV, med,
1781 logicalName.c_str(), 0, 0, 0);
1782 logicalV->SetVisAttributes(m_VisAttributes.back());
1783 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, posZ * CLHEP::cm), logicalV,
1784 physicalName.c_str(), m_logicalCDC, false, id);
1785
1786 }
1787
1788 void GeoCDCCreator::createTube2(const double rmin, const double rmax,
1789 const double phis, const double phie,
1790 const double thick, const double posZ,
1791 const int id, G4Material* med,
1792 const string& name)
1793 {
1794 const string solidName = "solid" + name;
1795 const string logicalName = "logical" + name;
1796 const string physicalName = "physical" + name;
1797 G4Tubs* solidV = new G4Tubs(solidName.c_str(),
1798 rmin * CLHEP::cm,
1799 rmax * CLHEP::cm,
1800 thick * CLHEP::cm / 2.0,
1801 phis,
1802 phie);
1803 G4LogicalVolume* logicalV = new G4LogicalVolume(solidV, med,
1804 logicalName.c_str(), 0, 0, 0);
1805 logicalV->SetVisAttributes(m_VisAttributes.back());
1806 new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, posZ * CLHEP::cm - thick * CLHEP::cm / 2.0), logicalV,
1807 physicalName.c_str(), m_logicalCDC, false, id);
1808
1809 }
1810
1811 void GeoCDCCreator::createMapper(G4LogicalVolume& topVolume)
1812 {
1814 if (!gcp.getMapperGeometry()) return;
1815
1816 const double xc = 0.5 * (-0.0002769 + 0.0370499) * CLHEP::cm;
1817 const double yc = 0.5 * (-0.0615404 + -0.108948) * CLHEP::cm;
1818 const double zc = 0.5 * (-35.3 + 48.5) * CLHEP::cm;
1819 //3 plates
1820 // const double plateWidth = 13.756 * CLHEP::cm;
1821 // const double plateThick = 1.203 * CLHEP::cm;
1822 // const double plateLength = 83.706 * CLHEP::cm;
1823 const double plateWidth = 13.8 * CLHEP::cm;
1824 const double plateThick = 1.2 * CLHEP::cm;
1825 const double plateLength = 83.8 * CLHEP::cm;
1826 const double phi = gcp.getMapperPhiAngle() * CLHEP::deg; //phi-angle in lab.
1827 // std::cout << "phi= " << phi << std::endl;
1828 // const double endRingRmin = 4.1135 * CLHEP::cm;
1829 // const double endRingRmax = 15.353 * CLHEP::cm;
1830 // const double endRingThick = 2.057 * CLHEP::cm;
1831 const double endPlateRmin = 4.0 * CLHEP::cm;
1832 const double endPlateRmax = 15.5 * CLHEP::cm;
1833 const double bwdEndPlateThick = 1.7 * CLHEP::cm;
1834 const double fwdEndPlateThick = 2.0 * CLHEP::cm;
1835
1836 G4Material* medAluminum = geometry::Materials::get("Al");
1837
1838 string name = "Plate";
1839 int pID = 0;
1840 G4Box* plateShape = new G4Box("solid" + name, .5 * plateWidth, .5 * plateThick, .5 * plateLength);
1841 G4LogicalVolume* logical0 = new G4LogicalVolume(plateShape, medAluminum, "logical" + name, 0, 0, 0);
1842 logical0->SetVisAttributes(m_VisAttributes.back());
1843 // const double x = .5 * plateWidth;
1844 const double x = xc + 0.5 * plateWidth;
1845 // const double y = endRingRmin;
1846 const double y = yc + endPlateRmin + 0.1 * CLHEP::cm;
1847 // double z = 2.871 * CLHEP::cm;
1848 G4ThreeVector xyz(x, y, zc);
1849 G4RotationMatrix rotM3 = G4RotationMatrix();
1850 xyz.rotateZ(phi);
1851 rotM3.rotateZ(phi);
1852 new G4PVPlacement(G4Transform3D(rotM3, xyz), logical0, "physical" + name, &topVolume, false, pID);
1853
1854 const double alf = 120. * CLHEP::deg;
1855 xyz.rotateZ(alf);
1856 rotM3.rotateZ(alf);
1857 new G4PVPlacement(G4Transform3D(rotM3, xyz), logical0, "physical" + name, &topVolume, false, pID + 1);
1858
1859 xyz.rotateZ(alf);
1860 rotM3.rotateZ(alf);
1861 new G4PVPlacement(G4Transform3D(rotM3, xyz), logical0, "physical" + name, &topVolume, false, pID + 2);
1862
1863 //Define 2 end-plates
1864 //bwd
1865 name = "BwdEndPlate";
1866 G4Tubs* BwdEndPlateShape = new G4Tubs("solid" + name, endPlateRmin, endPlateRmax, 0.5 * bwdEndPlateThick, 0., 360.*CLHEP::deg);
1867 G4LogicalVolume* logical1 = new G4LogicalVolume(BwdEndPlateShape, medAluminum, "logical" + name, 0, 0, 0);
1868 logical1->SetVisAttributes(m_VisAttributes.back());
1869 // z = -40.0105 * CLHEP::cm;
1870 double z = -35.3 * CLHEP::cm - 0.5 * bwdEndPlateThick;
1871 pID = 0;
1872 new G4PVPlacement(0, G4ThreeVector(xc, yc, z), logical1, "physical" + name, &topVolume, false, pID);
1873
1874 //fwd
1875 // z = 45.7525 * CLHEP::cm;
1876 // new G4PVPlacement(0, G4ThreeVector(0., 0., z), logical1, "physical" + name, &topVolume, false, pID + 1);
1877 name = "FwdEndPlate";
1878 G4Tubs* FwdEndPlateShape = new G4Tubs("solid" + name, endPlateRmin, endPlateRmax, 0.5 * fwdEndPlateThick, 0., 360.*CLHEP::deg);
1879 G4LogicalVolume* logical2 = new G4LogicalVolume(FwdEndPlateShape, medAluminum, "logical" + name, 0, 0, 0);
1880 logical2->SetVisAttributes(m_VisAttributes.back());
1881 z = 48.5 * CLHEP::cm + 0.5 * fwdEndPlateThick;
1882 new G4PVPlacement(0, G4ThreeVector(xc, yc, z), logical2, "physical" + name, &topVolume, false, pID);
1883 }
1884 }
1886}
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
void SetX(DataType x)
set X/1st-coordinate
Definition: B2Vector3.h:457
B2Vector3< DataType > Cross(const B2Vector3< DataType > &p) const
Cross product.
Definition: B2Vector3.h:296
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
void SetY(DataType y)
set Y/2nd-coordinate
Definition: B2Vector3.h:459
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:200
The Class for BeamBackground Sensitive Detector.
int getNEndPlateLayers() const
Get the number of endplate layers.
Definition: CDCGeometry.h:1382
double getR() const
Get Radius.
Definition: CDCGeometry.h:1112
double getZbwd() const
Get bwd z-position.
Definition: CDCGeometry.h:1122
std::vector< double > getRmin() const
Get the list of the Rmin corrdinates.
Definition: CDCGeometry.h:936
The Class for CDC geometry.
Definition: CDCGeometry.h:27
double getGlobalOffsetY() const
Get the global y offset of CDC wrt Belle2 coord.
Definition: CDCGeometry.h:1438
FieldLayer getFieldLayer(int i) const
Get the i-th field layer.
Definition: CDCGeometry.h:1553
OuterWall getOuterWall(int i) const
Get the i-th outer wall.
Definition: CDCGeometry.h:1483
std::vector< Rib2 > getRib2s() const
Get the list of rib2s.
Definition: CDCGeometry.h:1589
double getFiducialRmin() const
Get the fiducial Rmin of CDC sensitive volume.
Definition: CDCGeometry.h:1463
std::vector< Rib > getRibs() const
Get the list of ribs.
Definition: CDCGeometry.h:1584
double getGlobalOffsetX() const
Get the global x offset of CDC wrt Belle2 coord.
Definition: CDCGeometry.h:1433
std::vector< InnerWall > getInnerWalls() const
Get the list of inner walls.
Definition: CDCGeometry.h:1488
std::vector< Rib5 > getRib5s() const
Get the list of rib5s.
Definition: CDCGeometry.h:1604
double getFiducialRmax() const
Get the fiducial Rmax of CDC sensitive volume.
Definition: CDCGeometry.h:1468
int getNFieldWires() const
Get the number of field wires.
Definition: CDCGeometry.h:1615
std::vector< Rib3 > getRib3s() const
Get the list of rib3s.
Definition: CDCGeometry.h:1594
MotherVolume getMotherVolume() const
Get the mother volume geometry of CDC.
Definition: CDCGeometry.h:1473
int getNSenseWires() const
Get the number of sense wires.
Definition: CDCGeometry.h:1548
double getSenseDiameter() const
Get the diameter of sense wire.
Definition: CDCGeometry.h:1538
std::vector< Frontend > getFrontends() const
Get the list of frontend layers.
Definition: CDCGeometry.h:1518
InnerWall getInnerWall(int i) const
Get the i-th inner wall.
Definition: CDCGeometry.h:1493
double getFieldDiameter() const
Get the diameter of field wire.
Definition: CDCGeometry.h:1610
EndPlate getEndPlate(int i) const
Get the i-th endplate.
Definition: CDCGeometry.h:1503
std::vector< Rib4 > getRib4s() const
Get the list of rib4s.
Definition: CDCGeometry.h:1599
int getNSenseLayers() const
Get the number of sense layers.
Definition: CDCGeometry.h:1533
std::vector< EndPlate > getEndPlates() const
Get the list of endplates.
Definition: CDCGeometry.h:1508
double getFeedthroughLength() const
Get the length of feedthrough.
Definition: CDCGeometry.h:1620
double getGlobalOffsetZ() const
Get the global z offset of CDC wrt Belle2 coord.
Definition: CDCGeometry.h:1443
std::vector< OuterWall > getOuterWalls() const
Get the list of outer walls.
Definition: CDCGeometry.h:1478
SenseLayer getSenseLayer(int i) const
Get i-th sense layer.
Definition: CDCGeometry.h:1523
The Class for CDC Geometry Control Parameters.
bool getPrintMaterialTable() const
Get printMaterialTable flag.
double getMaterialDefinitionMode() const
Get material definition mode.
bool getMapperGeometry()
Get mapper geometry flag.
double getMapperPhiAngle()
Get mapper phi-angle.
static CDCGeoControlPar & getInstance()
Static method to get a reference to the CDCGeoControlPar instance.
The Class for CDC Geometry Parameters.
double fieldWireDiameter() const
Returns diameter of the field wire.
const B2Vector3D wireForwardPosition(uint layerId, int cellId, EWirePosition set=c_Base) const
Returns the forward position of the input sense wire.
double senseWireDiameter() const
Returns diameter of the sense wire.
const B2Vector3D wireBackwardPosition(uint layerId, int cellId, EWirePosition set=c_Base) const
Returns the backward position of the input sense wire.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
The Class for CDC Sensitive Detector.
static CDCSimControlPar & getInstance()
Static method to get a reference to the CDCSimControlPar instance.
void createBox(const double length, const double height, const double thick, const double x, const double y, const double z, const int id, G4Material *med, const std::string &name)
Create G4Box.
~GeoCDCCreator()
The destructor of the GeoCDCCreator class.
std::vector< G4VisAttributes * > m_VisAttributes
Vector of pointers to G4VisAttributes.
std::vector< G4UserLimits * > m_userLimits
Vector of pointers to G4UserLimits.
BkgSensitiveDetector * m_bkgsensitive
Sensitive detector for background studies.
void createGeometry(const CDCGeometry &parameters, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Create G4 geometry of CDC.
void createTube2(const double rmin, const double rmax, const double phis, const double phie, const double thick, const double posZ, const int id, G4Material *med, const std::string &name)
Create G4Tube2.
void createTube(const double rmin, const double rmax, const double thick, const double posZ, const int id, G4Material *med, const std::string &name)
Create G4Tube.
GeoCDCCreator()
Constructor of the GeoCDCCreator class.
void createNeutronShields(const GearDir &content)
Create neutron shield from gearbox.
G4VPhysicalVolume * m_physicalCDC
CDC G4 physical volume.
void createTorus(const double rmin1, const double rmax1, const double thick, const double posZ, const int id, G4Material *med, const std::string &name)
Create G4Torus.
std::vector< BkgSensitiveDetector * > m_BkgSensitiveRib4
Sensitive detectors for background studies (rib4).
void createCone(const double rmin1, const double rmax1, const double rmin2, const double rmax2, const double thick, const double posz, const int id, G4Material *med, const std::string &name)
Create G4Cone.
void createCover2s(const GearDir &content)
Create CDC cover2s from gear box.
void createMapper(G4LogicalVolume &topVolume)
Create the B-field mapper geometry (tentative function)
G4LogicalVolume * m_logicalCDC
CDC G4 logical volume.
void createCovers(const GearDir &content)
Create CDC covers from gear box.
CDCSensitiveDetector * m_sensitive
Sensitive detector.
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
void append(const std::string &path)
Append something to the current path, modifying the GearDir in place.
Definition: GearDir.h:52
virtual std::string getString(const std::string &path="") const noexcept(false) override
Get the parameter path as a string.
Definition: GearDir.h:69
double getLength(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard length unit.
Definition: Interface.h:259
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
geometry::CreatorFactory< GeoCDCCreator > GeoCDCFactory("CDCCreator")
Register the GeoCreator.
GeometryTypes
Flag indiciating the type of geometry to be used.
Abstract base class for different kinds of events.
STL namespace.
Very simple class to provide an easy way to register creators with the CreatorManager.