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