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