9 #include <ir/geometry/GeoCryostatCreator.h>
10 #include <ir/simulation/SensitiveDetector.h>
12 #include <geometry/Materials.h>
13 #include <geometry/CreatorFactory.h>
14 #include <geometry/utilities.h>
15 #include <framework/gearbox/Unit.h>
18 #include <boost/format.hpp>
19 #include <boost/foreach.hpp>
20 #include <boost/algorithm/string.hpp>
22 #include <G4LogicalVolume.hh>
23 #include <G4PVPlacement.hh>
27 #include <G4Polycone.hh>
28 #include <G4EllipticalTube.hh>
29 #include <G4UnionSolid.hh>
30 #include <G4IntersectionSolid.hh>
31 #include <G4SubtractionSolid.hh>
32 #include <G4UserLimits.hh>
35 using namespace boost;
44 using namespace geometry;
52 geometry::CreatorFactory<GeoCryostatCreator> GeoCryostatFactory(
"CryostatCreator");
58 GeoCryostatCreator::GeoCryostatCreator()
63 GeoCryostatCreator::~GeoCryostatCreator()
68 void GeoCryostatCreator::createGeometry(G4LogicalVolume& topVolume,
GeometryTypes)
147 double stepMax = 5.0 * Unit::mm;
148 int flag_limitStep = int(m_config.getParameter(
"LimitStepLength"));
150 const double unitFactor = Unit::cm / Unit::mm;
152 double crossingAngleHER = m_config.getParameter(
"CrossingAngle.HER", 0.0415);
153 double crossingAngleLER = m_config.getParameter(
"CrossingAngle.LER", -0.0415);
155 G4Transform3D transform_HER = G4Translate3D(0., 0., 0.);
156 transform_HER = transform_HER * G4RotateY3D(crossingAngleHER);
158 G4Transform3D transform_LER = G4Translate3D(0., 0., 0.);
159 transform_LER = transform_LER * G4RotateY3D(crossingAngleLER);
161 map<string, CryostatElement> elements;
168 std::string prep =
"TubeR.";
170 const int TubeR_N = int(m_config.getParameter(prep +
"N"));
172 std::vector<double> TubeR_Z(TubeR_N);
173 std::vector<double> TubeR_R(TubeR_N);
174 std::vector<double> TubeR_r(TubeR_N);
176 for (
int i = 0; i < TubeR_N; ++i) {
177 ostringstream ossZID;
180 ostringstream ossRID;
183 ostringstream ossrID;
186 TubeR_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
187 TubeR_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
188 TubeR_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
191 tubeR.
transform = G4Translate3D(0., 0., 0.);
192 tubeR.
geo =
new G4Polycone(
"geo_TubeR_name", 0, 2 * M_PI, TubeR_N, &(TubeR_Z[0]), &(TubeR_r[0]), &(TubeR_R[0]));
194 elements[
"TubeR"] = tubeR;
199 const int TubeR2_N = int(m_config.getParameter(prep +
"N"));
201 std::vector<double> TubeR2_Z(TubeR2_N);
202 std::vector<double> TubeR2_R(TubeR2_N);
203 std::vector<double> TubeR2_r(TubeR2_N);
205 for (
int i = 0; i < TubeR2_N; ++i) {
206 ostringstream ossZID;
209 ostringstream ossRID;
212 ostringstream ossrID;
215 TubeR2_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
216 TubeR2_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
217 TubeR2_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
220 tubeR2.
transform = G4Translate3D(0., 0., 0.);
221 tubeR2.
geo =
new G4Polycone(
"geo_TubeR2_name", 0, 2 * M_PI, TubeR2_N, &(TubeR2_Z[0]), &(TubeR2_r[0]), &(TubeR2_R[0]));
223 elements[
"TubeR2"] = tubeR2;
228 const int TubeL_N = int(m_config.getParameter(prep +
"N"));
230 std::vector<double> TubeL_Z(TubeL_N);
231 std::vector<double> TubeL_R(TubeL_N);
232 std::vector<double> TubeL_r(TubeL_N);
234 for (
int i = 0; i < TubeL_N; ++i) {
235 ostringstream ossZID;
238 ostringstream ossRID;
241 ostringstream ossrID;
244 TubeL_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
245 TubeL_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
246 TubeL_r[i] = m_config.getParameter(prep + ossrID.str()) * unitFactor;
249 tubeL.
transform = G4Translate3D(0., 0., 0.);
250 tubeL.
geo =
new G4Polycone(
"geo_TubeL_name", 0, 2 * M_PI, TubeL_N, &(TubeL_Z[0]), &(TubeL_r[0]), &(TubeL_R[0]));
252 elements[
"TubeL"] = tubeL;
263 const int A1spc1_N = int(m_config.getParameter(prep +
"N"));
265 std::vector<double> A1spc1_Z(A1spc1_N);
266 std::vector<double> A1spc1_r(A1spc1_N);
267 std::vector<double> A1spc1_R(A1spc1_N);
269 for (
int i = 0; i < A1spc1_N; ++i) {
270 ostringstream ossZID;
273 ostringstream ossRID;
276 ostringstream ossrID;
279 A1spc1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
280 A1spc1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
281 A1spc1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
285 G4Polycone* geo_A1spc1xx =
new G4Polycone(
"geo_A1spc1xx_name", 0, 2 * M_PI, A1spc1_N, &(A1spc1_Z[0]), &(A1spc1_r[0]),
291 const int A1spc2_N = int(m_config.getParameter(prep +
"N"));
293 std::vector<double> A1spc2_Z(A1spc2_N);
294 std::vector<double> A1spc2_R(A1spc2_N);
295 std::vector<double> A1spc2_r(A1spc2_N);
297 for (
int i = 0; i < A1spc2_N; ++i) {
298 ostringstream ossZID;
301 ostringstream ossRID;
304 ostringstream ossrID;
307 A1spc2_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
308 A1spc2_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
309 A1spc2_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
313 G4Polycone* geo_A1spc2xx =
new G4Polycone(
"geo_A1spc2xx_name", 0, 2 * M_PI, A1spc2_N, &(A1spc2_Z[0]), &(A1spc2_r[0]),
319 const int B1spc1_N = int(m_config.getParameter(prep +
"N"));
321 std::vector<double> B1spc1_Z(B1spc1_N);
322 std::vector<double> B1spc1_R(B1spc1_N);
323 std::vector<double> B1spc1_r(B1spc1_N);
325 for (
int i = 0; i < B1spc1_N; ++i) {
326 ostringstream ossZID;
329 ostringstream ossRID;
332 ostringstream ossrID;
335 B1spc1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
336 B1spc1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
337 B1spc1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
341 G4Polycone* geo_B1spc1xx =
new G4Polycone(
"geo_B1spc1xx_name", 0, 2 * M_PI, B1spc1_N, &(B1spc1_Z[0]), &(B1spc1_r[0]),
347 const int B1spc2_N = int(m_config.getParameter(prep +
"N"));
349 std::vector<double> B1spc2_Z(B1spc2_N);
350 std::vector<double> B1spc2_R(B1spc2_N);
351 std::vector<double> B1spc2_r(B1spc2_N);
353 for (
int i = 0; i < B1spc2_N; ++i) {
354 ostringstream ossZID;
357 ostringstream ossRID;
360 ostringstream ossrID;
363 B1spc2_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
364 B1spc2_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
365 B1spc2_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
369 G4Polycone* geo_B1spc2xx =
new G4Polycone(
"geo_B1spc2xx_name", 0, 2 * M_PI, B1spc2_N, &(B1spc2_Z[0]), &(B1spc2_r[0]),
373 B1spc2.
geo =
new G4IntersectionSolid(
"geo_B1spc2_name", geo_B1spc2xx, elements[
"TubeR2"].geo, B1spc2.
transform.inverse());
376 G4IntersectionSolid* geo_B1spc1x =
new G4IntersectionSolid(
"geo_B1spc1x_name", geo_B1spc1xx, elements[
"TubeR"].geo,
378 B1spc1.
geo =
new G4UnionSolid(
"geo_B1spc1_name", geo_B1spc1x, B1spc2.
geo);
380 A1spc2.
geo =
new G4IntersectionSolid(
"geo_A1spc2_name", geo_A1spc2xx, elements[
"TubeR2"].geo, A1spc2.
transform.inverse());
383 G4IntersectionSolid* geo_A1spc1xy =
new G4IntersectionSolid(
"geo_A1spc1xy_name", geo_A1spc1xx, elements[
"TubeR"].geo,
385 G4UnionSolid* geo_A1spc1x =
new G4UnionSolid(
"geo_A1spc1x_name", geo_A1spc1xy, A1spc2.
geo);
386 A1spc1.
geo =
new G4SubtractionSolid(
"geo_A1spc1_name", geo_A1spc1x, B1spc1.
geo, A1spc1.
transform.inverse()*B1spc1.
transform);
388 string strMat_A1spc1 = m_config.getParameterStr(
"A1spc1.Material");
389 G4Material* mat_A1spc1 = Materials::get(strMat_A1spc1);
390 A1spc1.
logi =
new G4LogicalVolume(A1spc1.
geo, mat_A1spc1,
"logi_A1spc1_name");
392 A1spc1.
logi->SetUserLimits(
new G4UserLimits(stepMax));
397 new G4PVPlacement(A1spc1.
transform, A1spc1.
logi,
"phys_A1spc1_name", &topVolume,
false, 0);
399 string strMat_B1spc1 = m_config.getParameterStr(
"B1spc1.Material");
400 G4Material* mat_B1spc1 = Materials::get(strMat_B1spc1);
401 B1spc1.
logi =
new G4LogicalVolume(B1spc1.
geo, mat_B1spc1,
"logi_B1spc1_name");
403 B1spc1.
logi->SetUserLimits(
new G4UserLimits(stepMax));
408 new G4PVPlacement(B1spc1.
transform, B1spc1.
logi,
"phys_B1spc1_name", &topVolume,
false, 0);
410 elements[
"A1spc1"] = A1spc1;
411 elements[
"A1spc2"] = A1spc2;
412 elements[
"B1spc1"] = B1spc1;
413 elements[
"B1spc2"] = B1spc2;
421 const int C1wal1_N = m_config.getParameter(prep +
"N");
423 std::vector<double> C1wal1_Z(C1wal1_N);
424 std::vector<double> C1wal1_R(C1wal1_N);
425 std::vector<double> C1wal1_r(C1wal1_N);
427 for (
int i = 0; i < C1wal1_N; ++i) {
428 ostringstream ossZID;
431 ostringstream ossRID;
434 ostringstream ossrID;
437 C1wal1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
438 C1wal1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
439 C1wal1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
442 C1wal1.
transform = G4Translate3D(0., 0., 0.);
445 G4Polycone* geo_C1wal1xxx =
new G4Polycone(
"geo_C1wal1xxx_name", 0, 2 * M_PI, C1wal1_N, &(C1wal1_Z[0]), &(C1wal1_r[0]),
447 G4IntersectionSolid* geo_C1wal1xx =
new G4IntersectionSolid(
"geo_C1wal1xx_name", geo_C1wal1xxx, elements[
"TubeR"].geo,
448 elements[
"TubeR"].transform);
449 G4SubtractionSolid* geo_C1wal1x =
new G4SubtractionSolid(
"geo_C1wal1x_name", geo_C1wal1xx, elements[
"A1spc1"].geo,
450 elements[
"A1spc1"].transform);
451 C1wal1.
geo =
new G4SubtractionSolid(
"geo_C1wal1_name", geo_C1wal1x, elements[
"B1spc1"].geo, elements[
"B1spc1"].transform);
453 string strMat_C1wal1 = m_config.getParameterStr(prep +
"Material");
454 G4Material* mat_C1wal1 = Materials::get(strMat_C1wal1);
455 C1wal1.
logi =
new G4LogicalVolume(C1wal1.
geo, mat_C1wal1,
"logi_C1wal1_name");
460 new G4PVPlacement(0, G4ThreeVector(0, 0, 0), C1wal1.
logi,
"phys_C1wal1_name", &topVolume,
false, 0);
462 elements[
"C1wal1"] = C1wal1;
470 const int D1spc1_N = m_config.getParameter(prep +
"N");
472 std::vector<double> D1spc1_Z(D1spc1_N);
473 std::vector<double> D1spc1_r(D1spc1_N);
474 std::vector<double> D1spc1_R(D1spc1_N);
476 for (
int i = 0; i < D1spc1_N; ++i) {
477 ostringstream ossZID;
480 ostringstream ossRID;
483 ostringstream ossrID;
486 D1spc1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
487 D1spc1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
488 D1spc1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
492 G4Polycone* geo_D1spc1xx =
new G4Polycone(
"geo_D1spc1xx_name", 0, 2 * M_PI, D1spc1_N, &(D1spc1_Z[0]), &(D1spc1_r[0]),
498 const int E1spc1_N = int(m_config.getParameter(prep +
"N"));
500 std::vector<double> E1spc1_Z(E1spc1_N);
501 std::vector<double> E1spc1_R(E1spc1_N);
502 std::vector<double> E1spc1_r(E1spc1_N);
504 for (
int i = 0; i < E1spc1_N; ++i) {
505 ostringstream ossZID;
508 ostringstream ossRID;
511 ostringstream ossrID;
514 E1spc1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
515 E1spc1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
516 E1spc1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
520 G4Polycone* geo_E1spc1xx =
new G4Polycone(
"geo_E1spc1xx_name", 0, 2 * M_PI, E1spc1_N, &(E1spc1_Z[0]), &(E1spc1_r[0]),
524 G4IntersectionSolid* geo_D1spc1x =
new G4IntersectionSolid(
"geo_D1spc1x_name", geo_D1spc1xx, elements[
"TubeL"].geo,
526 E1spc1.
geo =
new G4IntersectionSolid(
"geo_E1spc1_name", geo_E1spc1xx, elements[
"TubeL"].geo, E1spc1.
transform.inverse());
527 D1spc1.
geo =
new G4SubtractionSolid(
"geo_D1spc1_name", geo_D1spc1x, E1spc1.
geo, D1spc1.
transform.inverse()*E1spc1.
transform);
529 string strMat_D1spc1 = m_config.getParameterStr(
"D1spc1.Material");
530 G4Material* mat_D1spc1 = Materials::get(strMat_D1spc1);
531 D1spc1.
logi =
new G4LogicalVolume(D1spc1.
geo, mat_D1spc1,
"logi_D1spc1_name");
533 D1spc1.
logi->SetUserLimits(
new G4UserLimits(stepMax));
538 new G4PVPlacement(D1spc1.
transform, D1spc1.
logi,
"phys_D1spc1_name", &topVolume,
false, 0);
540 string strMat_E1spc1 = m_config.getParameterStr(prep +
"Material");
541 G4Material* mat_E1spc1 = Materials::get(strMat_E1spc1);
542 E1spc1.
logi =
new G4LogicalVolume(E1spc1.
geo, mat_E1spc1,
"logi_E1spc1_name");
544 E1spc1.
logi->SetUserLimits(
new G4UserLimits(stepMax));
549 new G4PVPlacement(E1spc1.
transform, E1spc1.
logi,
"phys_E1spc1_name", &topVolume,
false, 0);
551 elements[
"E1spc1"] = E1spc1;
552 elements[
"D1spc1"] = D1spc1;
561 const int F1wal1_N = int(m_config.getParameter(prep +
"N"));
563 std::vector<double> F1wal1_Z(F1wal1_N);
564 std::vector<double> F1wal1_R(F1wal1_N);
565 std::vector<double> F1wal1_r(F1wal1_N);
567 for (
int i = 0; i < F1wal1_N; ++i) {
568 ostringstream ossZID;
571 ostringstream ossRID;
574 ostringstream ossrID;
577 F1wal1_Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
578 F1wal1_R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
579 F1wal1_r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
582 F1wal1.
transform = G4Translate3D(0., 0., 0.);
585 G4Polycone* geo_F1wal1xxx =
new G4Polycone(
"geo_F1wal1xxx_name", 0, 2 * M_PI, F1wal1_N, &(F1wal1_Z[0]), &(F1wal1_r[0]),
587 G4IntersectionSolid* geo_F1wal1xx =
new G4IntersectionSolid(
"geo_F1wal1xx_name", geo_F1wal1xxx, elements[
"TubeL"].geo,
588 elements[
"TubeL"].transform);
589 G4SubtractionSolid* geo_F1wal1x =
new G4SubtractionSolid(
"geo_F1wal1x_name", geo_F1wal1xx, elements[
"D1spc1"].geo,
590 elements[
"D1spc1"].transform);
591 F1wal1.
geo =
new G4SubtractionSolid(
"geo_F1wal1_name", geo_F1wal1x, elements[
"E1spc1"].geo, elements[
"E1spc1"].transform);
593 string strMat_F1wal1 = m_config.getParameterStr(prep +
"Material");
594 G4Material* mat_F1wal1 = Materials::get(strMat_F1wal1);
595 F1wal1.
logi =
new G4LogicalVolume(F1wal1.
geo, mat_F1wal1,
"logi_F1wal1_name");
600 new G4PVPlacement(F1wal1.
transform, F1wal1.
logi,
"phys_F1wal1_name", &topVolume,
false, 0);
602 elements[
"F1wal1"] = F1wal1;
607 std::vector<std::string> straightSections;
608 boost::split(straightSections, m_config.getParameterStr(
"Straight"), boost::is_any_of(
" "));
609 for (
const auto& name : straightSections) {
614 int N = int(m_config.getParameter(prep +
"N"));
616 std::vector<double> Z(N);
617 std::vector<double> R(N);
618 std::vector<double> r(N);
620 for (
int i = 0; i < N; ++i) {
621 ostringstream ossZID;
624 ostringstream ossRID;
627 ostringstream ossrID;
630 Z[i] = m_config.getParameter(prep + ossZID.str()) * unitFactor;
631 R[i] = m_config.getParameter(prep + ossRID.str()) * unitFactor;
632 r[i] = m_config.getParameter(prep + ossrID.str(), 0.0) * unitFactor;
635 polycone.
transform = G4Translate3D(0.0, 0.0, 0.0);
638 string motherVolume = m_config.getParameterStr(prep +
"MotherVolume");
639 string subtract = m_config.getParameterStr(prep +
"Subtract",
"");
640 string intersect = m_config.getParameterStr(prep +
"Intersect",
"");
642 string geo_polyconexx_name =
"geo_" + name +
"xx_name";
643 string geo_polyconex_name =
"geo_" + name +
"x_name";
644 string geo_polycone_name =
"geo_" + name +
"_name";
646 G4VSolid* geo_polyconexx, *geo_polycone;
649 geo_polyconexx =
new G4Polycone(geo_polyconexx_name, 0.0, 2 * M_PI, N, &(Z[0]), &(r[0]), &(R[0]));
650 G4VSolid* geo_polyconex =
new G4SubtractionSolid(geo_polyconex_name, geo_polyconexx, elements[subtract].geo,
651 elements[motherVolume].transform.inverse()*polycone.
transform.inverse()*elements[subtract].transform);
652 geo_polycone =
new G4IntersectionSolid(geo_polycone_name, geo_polyconex, elements[
intersect].geo,
653 elements[motherVolume].transform.inverse()*polycone.
transform.inverse()*elements[
intersect].transform);
654 }
else if (subtract !=
"") {
655 geo_polyconexx =
new G4Polycone(geo_polyconexx_name, 0.0, 2 * M_PI, N, &(Z[0]), &(r[0]), &(R[0]));
656 geo_polycone =
new G4SubtractionSolid(geo_polycone_name, geo_polyconexx, elements[subtract].geo,
657 elements[motherVolume].transform.inverse()*polycone.
transform.inverse()*elements[subtract].transform);
659 geo_polyconexx =
new G4Polycone(geo_polyconexx_name, 0.0, 2 * M_PI, N, &(Z[0]), &(r[0]), &(R[0]));
660 geo_polycone =
new G4IntersectionSolid(geo_polycone_name, geo_polyconexx, elements[
intersect].geo,
661 elements[motherVolume].transform.inverse()*polycone.
transform.inverse()*elements[
intersect].transform);
663 geo_polycone =
new G4Polycone(geo_polycone_name, 0.0, 2 * M_PI, N, &(Z[0]), &(r[0]), &(R[0]));
665 polycone.
geo = geo_polycone;
668 string strMat_polycone = m_config.getParameterStr(prep +
"Material");
669 G4Material* mat_polycone = Materials::get(strMat_polycone);
670 string logi_polycone_name =
"logi_" + name +
"_name";
671 polycone.
logi =
new G4LogicalVolume(polycone.
geo, mat_polycone, logi_polycone_name);
676 string phys_polycone_name =
"phys_" + name +
"_name";
677 new G4PVPlacement(polycone.
transform, polycone.
logi, phys_polycone_name, elements[motherVolume].logi,
false, 0);
682 elements[name] = polycone;
686 G4Tubs* geo_rvcR =
new G4Tubs(
"geo_rvcR", 60, 60 + 60, (620 - 560) / 2., 0, 2 * M_PI);
687 G4LogicalVolume* logi_rvcR =
new G4LogicalVolume(geo_rvcR, Materials::get(
"SUS316L"),
"logi_rvcR_name");
688 new G4PVPlacement(0, G4ThreeVector(0, 0, (620 + 560) / 2.), logi_rvcR,
"phys_rvcR_name", &topVolume,
false, 0);
690 G4Tubs* geo_rvcL =
new G4Tubs(
"geo_rvcL", 60, 60 + 60, (-560 - (-620)) / 2., 0, 2 * M_PI);
691 G4LogicalVolume* logi_rvcL =
new G4LogicalVolume(geo_rvcL, Materials::get(
"SUS316L"),
"logi_rvcL_name");
692 new G4PVPlacement(0, G4ThreeVector(0, 0, (-620 - 560) / 2.), logi_rvcL,
"phys_rvcL_name", &topVolume,
false, 0);
696 G4EllipticalTube* geo_elp_QC1LEx =
new G4EllipticalTube(
"geo_elp_QC1LEx", 10.5, 13.5, (-675 - (-1225)) / 2.);
697 G4IntersectionSolid* geo_elp_QC1LE =
new G4IntersectionSolid(
"geo_elp_QC1LE", elements[
"D2wal1"].geo, geo_elp_QC1LEx,
698 G4Translate3D(0, 0, (-675 - 1225) / 2.));
699 G4LogicalVolume* logi_elp_QC1LE =
new G4LogicalVolume(geo_elp_QC1LE, Materials::get(
"Vacuum"),
"logi_elp_QC1LE_name");
700 new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logi_elp_QC1LE,
"phys_elp_QC1LE_name", elements[
"D2wal1"].logi,
false, 0);
702 G4EllipticalTube* geo_elp_QC1LPx =
new G4EllipticalTube(
"geo_elp_QC1LPx", 10.5, 13.5, (-675 - (-1225)) / 2.);
703 G4IntersectionSolid* geo_elp_QC1LP =
new G4IntersectionSolid(
"geo_elp_QC1LP", elements[
"E2wal1"].geo, geo_elp_QC1LPx,
704 G4Translate3D(0, 0, (-675 - 1225) / 2.));
705 G4LogicalVolume* logi_elp_QC1LP =
new G4LogicalVolume(geo_elp_QC1LP, Materials::get(
"Vacuum"),
"logi_elp_QC1LP_name");
706 new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logi_elp_QC1LP,
"phys_elp_QC1LP_name", elements[
"E2wal1"].logi,
false, 0);
708 G4EllipticalTube* geo_elp_QC1REx =
new G4EllipticalTube(
"geo_elp_QC1REx", 10.5, 13.5, (1225 - 675) / 2.);
709 G4IntersectionSolid* geo_elp_QC1RE =
new G4IntersectionSolid(
"geo_elp_QC1RE", elements[
"A2wal1"].geo, geo_elp_QC1REx,
710 G4Translate3D(0, 0, (1225 + 675) / 2.));
711 G4LogicalVolume* logi_elp_QC1RE =
new G4LogicalVolume(geo_elp_QC1RE, Materials::get(
"Vacuum"),
"logi_elp_QC1RE_name");
712 new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logi_elp_QC1RE,
"phys_elp_QC1RE_name", elements[
"A2wal1"].logi,
false, 0);
714 G4EllipticalTube* geo_elp_QC1RPx =
new G4EllipticalTube(
"geo_elp_QC1RPx", 10.5, 13.5, (1225 - 675) / 2.);
715 G4IntersectionSolid* geo_elp_QC1RP =
new G4IntersectionSolid(
"geo_elp_QC1RP", elements[
"B2wal1"].geo, geo_elp_QC1RPx,
716 G4Translate3D(0, 0, (1225 + 675) / 2.));
717 G4LogicalVolume* logi_elp_QC1RP =
new G4LogicalVolume(geo_elp_QC1RP, Materials::get(
"Vacuum"),
"logi_elp_QC1RP_name");
718 new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logi_elp_QC1RP,
"phys_elp_QC1RP_name", elements[
"B2wal1"].logi,
false, 0);
The IR Sensitive Detector class.
int intersect(const TRGCDCLpar &lp1, const TRGCDCLpar &lp2, CLHEP::HepVector &v1, CLHEP::HepVector &v2)
intersection
void setVisibility(G4LogicalVolume &volume, bool visible)
Helper function to quickly set the visibility of a given volume.
void setColor(G4LogicalVolume &volume, const std::string &color)
Set the color of a logical volume.
GeometryTypes
Flag indiciating the type of geometry to be used.
Abstract base class for different kinds of events.
The struct for CryostatElement.
G4LogicalVolume * logi
Logical volume.
G4VSolid * geo
Solid volume.
G4Transform3D transform
Transformation.