Belle II Software development
GeoECLCreator Class Reference

The GeoECLCreator class. More...

#include <GeoECLCreator.h>

Inheritance diagram for GeoECLCreator:
CreatorBase

Public Member Functions

 GeoECLCreator ()
 Constructor of the GeoECLCreator class.
 
 ~GeoECLCreator ()
 The destructor of the GeoECLCreator class.
 
virtual void create (const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type) override
 Function to actually create the geometry, has to be overridden by derived classes.
 
virtual void createFromDB (const std::string &name, G4LogicalVolume &topVolume, geometry::GeometryTypes type) override
 Function to create the geometry from the Database.
 
virtual void createPayloads (const GearDir &content, const IntervalOfValidity &iov) override
 Function to create the geometry database.
 
 BELLE2_DEFINE_EXCEPTION (DBNotImplemented, "Cannot create geometry from Database.")
 Exception that will be thrown in createFromDB if member is not yet implemented by creator.
 

Private Member Functions

void barrel (G4LogicalVolume &)
 Make the ECL barrel and then place elements inside it.
 
void backward (G4LogicalVolume &)
 Place elements inside the backward endcap.
 
void forward (G4LogicalVolume &)
 Place elements inside the forward endcap.
 
G4LogicalVolume * wrapped_crystal (const shape_t *s, const std::string &endcap, double wrapthickness)
 Wrapped crystal.
 
void defineVisAttributes ()
 Define visual attributes.
 
const G4VisAttributes * att (const std::string &n) const
 Define visual attributes.
 
G4LogicalVolume * get_preamp () const
 Get Logical volume of preamplifier.
 
double get_pa_box_height () const
 Getter for preamplifier box height (hard-coded to be 2)
 

Private Attributes

const ECLCrystalsShapeAndPositionm_sap
 pointer to a storage with crystal shapes and positions
 
Simulation::SensitiveDetectorBasem_sensitive
 Sensitive detector.
 
Simulation::SensitiveDetectorBasem_sensediode
 Sensitive diode.
 
std::map< std::string, G4VisAttributes * > m_atts
 Vector of background-Sensitive detectors.
 
int m_overlap
 overlap
 

Detailed Description

The GeoECLCreator class.

The creator for the ECL geometry of the Belle II detector.

Definition at line 34 of file GeoECLCreator.h.

Constructor & Destructor Documentation

◆ GeoECLCreator()

Constructor of the GeoECLCreator class.

Definition at line 43 of file GeoECLCreator.cc.

43 : m_sap(0), m_overlap(0)
44{
45 m_sensediode = nullptr;
46 m_sensitive = new SensitiveDetector("ECLSensitiveDetector", (2 * 24)*CLHEP::eV, 10 * CLHEP::MeV);
47 G4SDManager::GetSDMpointer()->AddNewDetector(m_sensitive);
49}
void defineVisAttributes()
Define visual attributes.
Simulation::SensitiveDetectorBase * m_sensediode
Sensitive diode.
Definition: GeoECLCreator.h:96
const ECLCrystalsShapeAndPosition * m_sap
pointer to a storage with crystal shapes and positions
Definition: GeoECLCreator.h:91
Simulation::SensitiveDetectorBase * m_sensitive
Sensitive detector.
Definition: GeoECLCreator.h:94
The Class for ARICH Sensitive Detector.

◆ ~GeoECLCreator()

The destructor of the GeoECLCreator class.

Definition at line 52 of file GeoECLCreator.cc.

53{
54 for (auto a : m_atts) delete a.second;
55}
std::map< std::string, G4VisAttributes * > m_atts
Vector of background-Sensitive detectors.
Definition: GeoECLCreator.h:98

Member Function Documentation

◆ att()

const G4VisAttributes * att ( const std::string &  n) const
private

Define visual attributes.

Parameters
nAttribute name

Definition at line 141 of file GeoECLCreator.cc.

142{
143 auto p = m_atts.find(n);
144 assert(p != m_atts.end());
145 return p->second;
146}

◆ backward()

void backward ( G4LogicalVolume &  _top)
private

Place elements inside the backward endcap.

Definition at line 42 of file backward.cc.

43{
44 G4LogicalVolume* top = &_top;
45
46 const bool sec = 0;
47 const double phi0 = 0;
48 const double dphi = sec ? M_PI / 16 : 2 * M_PI;
49
50 const bool b_inner_support_ring = true;
51 const bool b_outer_support_ring = true;
52 const bool b_support_wall = true;
53 const bool b_ribs = true;
54 const bool b_septum_wall = true;
55 const bool b_crystals = true;
56 const bool b_preamplifier = true;
57 const bool b_support_leg = true;
58 const int overlap = m_overlap;
59
60 int npoints = 1000 * 1000;
61
62 // vector<cplacement_t> bp = load_placements("/ecl/data/crystal_placement_backward.dat");
63 vector<cplacement_t> bp = load_placements(m_sap, ECLParts::backward);
64 vector<cplacement_t>::iterator fp = find_if(bp.begin(), bp.end(), [](const cplacement_t& p) {
65 const int ECL_backward_part = 1000; // lookup part
66 return p.nshape == ECL_backward_part;
67 });
68 // global transformation before placing the whole backward part in the top logical volume
69 G4Transform3D gTrans = (fp == bp.end()) ? G4Translate3D(0, 0, -1020) : get_transform(*fp);
70 if (fp != bp.end()) bp.erase(fp); // now not needed
71
72 if (b_inner_support_ring) {
73 zr_t vc1[] = {{0., 452.3 + 3}, {0., 452.3}, {3., 474.9 - 20 / cosd(27.81)}, {434., 702.27 - 20 / cosd(27.81)}, {434., 702.27}, {3., 474.9}, {3., 452.3 + 3}};
74 std::vector<zr_t> contour1(vc1, vc1 + sizeof(vc1) / sizeof(zr_t));
75 G4VSolid* part1solid = new BelleLathe("part1solid", phi0, dphi, contour1);
76 G4LogicalVolume* part1logical = new G4LogicalVolume(part1solid, Materials::get("SUS304"), "part1logical", 0, 0, 0);
77 part1logical->SetVisAttributes(att("iron"));
78 auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part1logical, "ECL_BackwardSupport_part1physical", top, false, 0, 0);
79 if (overlap) pv->CheckOverlaps(npoints);
80 }
81
82 if (b_support_wall) {
83 // solving equation to get L : 3+L*cosd(52.90)+1.6*cosd(52.90+90) = 435 - 202.67
84 double L = (435 - 202.67 - 3 - 1.6 * cosd(52.90 + 90)) / cosd(52.90);
85 zr_t vc23[] = {{0, 452.3 + 3}, {3, 452.3 + 3}, {3, 1190.2}, {3 + L * cosd(52.90), 1190.2 + L * sind(52.90)},
86 {3 + L * cosd(52.90) + 1.6 * cosd(52.90 + 90), 1190.2 + L * sind(52.90) + 1.6 * sind(52.90 + 90)}, {3 + 1.6 * cosd(52.90 + 90), 1190.2 + 1.6 * sind(52.90 + 90)}, {0, 1190.2}
87 };
88 std::vector<zr_t> contour23(vc23, vc23 + sizeof(vc23) / sizeof(zr_t));
89 G4VSolid* part23solid = new BelleLathe("part23solid", phi0, dphi, contour23);
90 G4LogicalVolume* part23logical = new G4LogicalVolume(part23solid, Materials::get("A5052"), "part23logical", 0, 0, 0);
91 part23logical->SetVisAttributes(att("alum"));
92 auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part23logical, "ECL_BackwardSupport_part23physical", top, false, 0, 0);
93 if (overlap) pv->CheckOverlaps(npoints);
94 }
95
96 if (b_outer_support_ring) {
97 zr_t vc4[] = {{434 - 214.8, 1496 - 20}, {434, 1496 - 20}, {434, 1496 - 5}, {434 + 5, 1496 - 5}, {434 + 5, 1496}, {434 - 199.66, 1496}};
98 std::vector<zr_t> contour4(vc4, vc4 + sizeof(vc4) / sizeof(zr_t));
99 G4VSolid* part4solid = new BelleLathe("part4solid", phi0, dphi, contour4);
100 G4LogicalVolume* part4logical = new G4LogicalVolume(part4solid, Materials::get("SUS304"), "part4logical", 0, 0, 0);
101 part4logical->SetVisAttributes(att("iron"));
102 auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part4logical, "ECL_BackwardSupport_part4physical", top, false, 0, 0);
103 if (overlap) pv->CheckOverlaps(npoints);
104 }
105
106 zr_t cont_array_in[] = {{3., 474.9}, {434., 702.27}, {434, 1496 - 20}, {434 - 214.8, 1496 - 20}, {3, 1190.2}};
107 std::vector<zr_t> contour_in(cont_array_in, cont_array_in + sizeof(cont_array_in) / sizeof(zr_t));
108 G4VSolid* innervolume_solid = new BelleLathe("innervolume_solid", 0, 2 * M_PI, contour_in);
109 G4LogicalVolume* innervolume_logical = new G4LogicalVolume(innervolume_solid, Materials::get("G4_AIR"),
110 "innervolume_logical", 0, 0, 0);
111 innervolume_logical->SetVisAttributes(att("air"));
112
113 // Set up region for production cuts
114 G4Region* aRegion = new G4Region("ECLBackwardEnvelope");
115 innervolume_logical->SetRegion(aRegion);
116 aRegion->AddRootLogicalVolume(innervolume_logical);
117
118 auto gpvbp = new G4PVPlacement(gTrans * G4RotateY3D(M_PI),
119 innervolume_logical, "ECLBackwardPhysical",
120 top, false, 0, 0);
121 if (overlap) gpvbp->CheckOverlaps(npoints);
122
123 G4VSolid* innervolumesector_solid = new BelleLathe("innervolumesector_solid", -M_PI / 8, M_PI / 4, contour_in);
124 G4LogicalVolume* innervolumesector_logical = new G4LogicalVolume(innervolumesector_solid, Materials::get("G4_AIR"),
125 "innervolumesector_logical", 0, 0, 0);
126 innervolumesector_logical->SetVisAttributes(att("air"));
127 new G4PVReplica("ECLBackwardSectorPhysical", innervolumesector_logical, innervolume_logical, kPhi, 8, M_PI / 4, 0);
128
129 if (b_ribs) {
130 double H = 60, W = 20;
131 double X0 = 702.27 + 0.001, X1 = 1496 - 20;
132 G4TwoVector r0o(X1, 0), r1o(X1 * sqrt(1 - pow(W / X1, 2)), W);
133 double beta = asin(W / X0);
134 G4TwoVector r0i(X0 / cos(beta / 2), 0), r1i(X0 * cos(beta / 2) - tan(beta / 2) * (W - X0 * sin(beta / 2)), W);
135 double dxymzp = (r0o - r0i).x(), dxypzp = (r1o - r1i).x();
136 double theta = atan(tand(27.81) / 2);
137 double dxymzm = dxymzp + tand(27.81) * H, dxypzm = dxypzp + tand(27.81) * H;
138
139 G4TwoVector m0 = (r0i + r0o) * 0.5, m1 = (r1i + r1o) * 0.5, dm = m1 - m0;
140 double alpha = atan(dm.x() / dm.y());
141
142 G4VSolid* solid6_p1 = new G4Trap("solid6_p1", H / 2, theta, 0, W / 2, dxymzm / 2, dxypzm / 2, alpha, W / 2, dxymzp / 2, dxypzp / 2,
143 alpha);
144 G4LogicalVolume* lsolid6_p1 = new G4LogicalVolume(solid6_p1, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
145 lsolid6_p1->SetVisAttributes(att("asolid6"));
146 G4Transform3D tsolid6_p1(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, W / 2, 434 - H / 2));
147 auto pv61 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*tsolid6_p1, lsolid6_p1, "psolid6_p1", innervolumesector_logical, false, 0, 0);
148 if (overlap) pv61->CheckOverlaps(npoints);
149 auto pv62 = new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p1, lsolid6_p1, "psolid6_p2", innervolumesector_logical, false, 0, 0);
150 if (overlap) pv62->CheckOverlaps(npoints);
151
152 H = 40;
153 dxymzm = dxymzp + tand(27.81) * H, dxypzm = dxypzp + tand(27.81) * H;
154 G4VSolid* solid6_p2 = new G4Trap("solid6_p2", H / 2, theta, 0, W / 2, dxypzm / 2, dxymzm / 2, -alpha, W / 2, dxypzp / 2, dxymzp / 2,
155 -alpha);
156 G4LogicalVolume* lsolid6_p2 = new G4LogicalVolume(solid6_p2, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
157 lsolid6_p2->SetVisAttributes(att("asolid6"));
158 G4Transform3D tsolid6_p2(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, -W / 2, 434 - H / 2));
159 auto pv63 = new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p2, lsolid6_p2, "psolid6_p3", innervolumesector_logical, false, 0, 0);
160 if (overlap) pv63->CheckOverlaps(npoints);
161 auto pv64 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid6_p2, lsolid6_p2, "psolid6_p4", innervolumesector_logical, false, 0, 0);
162 if (overlap) pv64->CheckOverlaps(npoints);
163
164 G4VSolid* solid7_p8 = new G4Box("solid7_p8", 171. / 2, (140. - 40) / 2 / 2, 40. / 2);
165 G4LogicalVolume* lsolid7 = new G4LogicalVolume(solid7_p8, Materials::get("SUS304"), "lsolid7", 0, 0, 0);
166 lsolid7->SetVisAttributes(att("asolid7"));
167 double dx = sqrt(X1 * X1 - 70 * 70) - 171. / 2;
168 G4Transform3D tsolid7_p1(G4Translate3D(dx, -20 - 25, 434 - 40. / 2));
169 auto pv71 = new G4PVPlacement(tsolid7_p1, lsolid7, "psolid7_p1", innervolumesector_logical, false, 0, 0);
170 if (overlap) pv71->CheckOverlaps(npoints);
171 G4Transform3D tsolid7_p2(G4Translate3D(dx, 20 + 25, 434 - 40. / 2));
172 auto pv72 = new G4PVPlacement(tsolid7_p2, lsolid7, "psolid7_p2", innervolumesector_logical, false, 0, 0);
173 if (overlap) pv72->CheckOverlaps(npoints);
174
175 double L = X1 - (X0 - tand(27.81) * 40) - 10;
176 G4VSolid* solid13 = new G4Box("solid13", L / 2, 5. / 2, 18. / 2);
177 G4LogicalVolume* lsolid13 = new G4LogicalVolume(solid13, Materials::get("SUS304"), "lsolid13", 0, 0, 0);
178 lsolid13->SetVisAttributes(att("asolid13"));
179 G4Transform3D tsolid13(G4TranslateZ3D(434 - 60 + 18. / 2)*G4TranslateY3D(-5. / 2 - 0.5 / 2)*G4TranslateX3D(X0 - tand(
180 27.81) * 40 + L / 2 + 5));
181 auto pv131 = new G4PVPlacement(tsolid13, lsolid13, "psolid13_p1", innervolumesector_logical, false, 0, 0);
182 if (overlap) pv131->CheckOverlaps(npoints);
183 auto pv132 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid13, lsolid13, "psolid13_p2", innervolumesector_logical, false, 0, 0);
184 if (overlap) pv132->CheckOverlaps(npoints);
185 }
186
187
188 double zsep = 135;
189 if (b_septum_wall) {
190 double d = 5;
191 Point_t vin[] = {{434. - zsep, 702.27 - tand(27.81)* zsep}, {434. - 60, 702.27 - tand(27.81) * 60}, {434. - 60, 1496 - 20 - d}, {434. - zsep, 1496 - 20 - d}};
192 const int n = sizeof(vin) / sizeof(Point_t);
193 Point_t c = centerofgravity(vin, vin + n);
194 G4ThreeVector contour_swall[n * 2];
195 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2);
196 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2);
197
198 G4VSolid* septumwall_solid = new BelleCrystal("septumwall_solid", n, contour_swall);
199
200 G4LogicalVolume* septumwall_logical = new G4LogicalVolume(septumwall_solid, Materials::get("A5052"),
201 "septumwall_logical", 0, 0, 0);
202 septumwall_logical->SetVisAttributes(att("alum2"));
203 auto pv = new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0), septumwall_logical,
204 "septumwall_physical", innervolumesector_logical, false, 0, 0);
205 if (overlap) pv->CheckOverlaps(npoints);
206
207 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
208 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
209
210 G4VSolid* septumwall2_solid = new BelleCrystal("septumwall2_solid", n, contour_swall);
211
212 G4LogicalVolume* septumwall2_logical = new G4LogicalVolume(septumwall2_solid, Materials::get("A5052"),
213 "septumwall2_logical", 0, 0, 0);
214 septumwall2_logical->SetVisAttributes(att("alum2"));
215 auto pv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
216 0.5 / 2 / 2),
217 septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 0, 0);
218 if (overlap) pv0->CheckOverlaps(npoints);
219 auto pv1 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
220 -0.5 / 2 / 2),
221 septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 1, 0);
222 if (overlap) pv1->CheckOverlaps(npoints);
223 }
224
225 zr_t vcr[] = {{3., 474.9}, {434. - zsep, 702.27 - tand(27.81)* zsep}, {434 - zsep, 1496 - 20}, {434 - 214.8, 1496 - 20}, {3, 1190.2}};
226 std::vector<zr_t> ccr(vcr, vcr + sizeof(vcr) / sizeof(zr_t));
227 G4VSolid* crystalvolume_solid = new BelleLathe("crystalvolume_solid", 0, M_PI / 8, ccr);
228 G4LogicalVolume* crystalvolume_logical = new G4LogicalVolume(crystalvolume_solid, Materials::get("G4_AIR"),
229 "crystalvolume_logical", 0, 0, 0);
230 crystalvolume_logical->SetVisAttributes(att("air"));
231 auto gpv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8), crystalvolume_logical, "ECLBackwardCrystalSectorPhysical_0",
232 innervolumesector_logical,
233 false, 0, 0);
234 if (overlap) gpv0->CheckOverlaps(npoints);
235 auto gpv1 = new G4PVPlacement(G4RotateZ3D(0), crystalvolume_logical, "ECLBackwardCrystalSectorPhysical_1",
236 innervolumesector_logical, false, 1,
237 0);
238 if (overlap) gpv1->CheckOverlaps(npoints);
239
240 if (b_septum_wall) {
241 double d = 5, dr = 0.001;
242 Point_t vin[] = {{3., 474.9}, {434. - zsep, 702.27 - tand(27.81)* zsep}, {434 - zsep, 1496 - 20 - d - dr}, {434 - 214.8 - d / tand(52.90), 1496 - 20 - d - dr}, {3, 1190.2 - dr}};
243 const int n = sizeof(vin) / sizeof(Point_t);
244 Point_t c = centerofgravity(vin, vin + n);
245 G4ThreeVector contour_swall[n * 2];
246
247 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
248 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
249
250 G4VSolid* septumwall3_solid = new BelleCrystal("septumwall3_solid", n, contour_swall);
251
252 G4LogicalVolume* septumwall3_logical = new G4LogicalVolume(septumwall3_solid, Materials::get("A5052"),
253 "septumwall3_logical", 0, 0, 0);
254 septumwall3_logical->SetVisAttributes(att("alum2"));
255 auto pv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0.5 / 2 / 2),
256 septumwall3_logical,
257 "septumwall3_physical_0", crystalvolume_logical, false, 0, overlap);
258 if (overlap) pv0->CheckOverlaps(npoints);
259 auto pv1 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
260 -0.5 / 2 / 2),
261 septumwall3_logical, "septumwall3_physical_1", crystalvolume_logical, false, 1, overlap);
262 if (overlap) pv1->CheckOverlaps(npoints);
263 }
264
265 if (b_crystals) {
266 // vector<shape_t*> cryst = load_shapes("/ecl/data/crystal_shape_backward.dat");
267 vector<shape_t*> cryst = load_shapes(m_sap, ECLParts::backward);
268 vector<G4LogicalVolume*> wrapped_crystals;
269 for (auto it = cryst.begin(); it != cryst.end(); it++) {
270 shape_t* s = *it;
271 wrapped_crystals.push_back(wrapped_crystal(s, "backward", 0.20 - 0.02));
272 }
273 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
274 const cplacement_t& t = *it;
275 auto s = find_if(cryst.begin(), cryst.end(), [&t](const shape_t* shape) {return shape->nshape == t.nshape;});
276 if (s == cryst.end()) continue;
277
278 G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(t);
279 int indx = it - bp.begin();
280 auto pv = new G4PVPlacement(twc, wrapped_crystals[s - cryst.begin()], suf("ECLBackwardWrappedCrystal_Physical", indx),
281 crystalvolume_logical,
282 false, (ECLElementNumbers::c_NCrystalsForwardBarrel) / 16 + indx, 0);
283 if (overlap)pv->CheckOverlaps(npoints);
284 }
285
286 for (shape_t* shape : cryst) {
287 delete shape;
288 }
289 }
290
291 if (b_preamplifier) {
292 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
293 G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(*it);
294 int indx = it - bp.begin();
295 auto pv = new G4PVPlacement(twc * G4TranslateZ3D(300 / 2 + 0.20 + get_pa_box_height() / 2)*G4RotateZ3D(-M_PI / 2), get_preamp(),
296 suf("phys_backward_preamplifier", indx), crystalvolume_logical, false, (ECLElementNumbers::c_NCrystalsForwardBarrel) / 16 + indx,
297 0);
298 if (overlap)pv->CheckOverlaps(npoints);
299 }
300 }
301
302 if (b_support_leg) {
303 const G4VisAttributes* batt = att("iron");
304
305 G4VSolid* s1 = new G4Box("leg_p1", 130. / 2, 185. / 2, (40. - 5) / 2);
306 G4LogicalVolume* l1 = new G4LogicalVolume(s1, Materials::get("SUS304"), "l1", 0, 0, 0);
307 G4Transform3D t1 = G4Translate3D(0, 185. / 2, (40. - 5) / 2);
308 l1->SetVisAttributes(batt);
309
310 Point_t v3[] = {{ -212. / 2, -135. / 2}, {212. / 2 - 30, -135. / 2}, {212. / 2, -135. / 2 + 30}, {212. / 2, 135. / 2}, { -212. / 2, 135. / 2}};
311 const int n3 = sizeof(v3) / sizeof(Point_t);
312 G4ThreeVector c3[n3 * 2];
313
314 for (int i = 0; i < n3; i++) c3[i + 0] = G4ThreeVector(v3[i].x, v3[i].y, -60. / 2);
315 for (int i = 0; i < n3; i++) c3[i + n3] = G4ThreeVector(v3[i].x, v3[i].y, 60. / 2);
316
317 G4VSolid* s3 = new BelleCrystal("leg_p3", n3, c3);
318 G4LogicalVolume* l3 = new G4LogicalVolume(s3, Materials::get("SUS304"), "l3", 0, 0, 0);
319 G4Transform3D t3 = G4Translate3D(0, 135. / 2 + 35, 40. - 5. + 212. / 2) * G4RotateY3D(-M_PI / 2);
320 l3->SetVisAttributes(batt);
321
322 G4VSolid* s6 = new G4Box("leg_p6", 140. / 2, 189. / 2, 160. / 2);
323 G4LogicalVolume* l6 = new G4LogicalVolume(s6, Materials::get("G4_AIR"), "l6", 0, 0, 0);
324 G4Transform3D t6 = G4Translate3D(0, 170. + 189. / 2, 57. + 35. + 160. / 2);
325 l6->SetVisAttributes(att("air"));
326
327 G4VSolid* s6a = new G4Box("leg_p6a", 140. / 2, (189. - 45.) / 2, 160. / 2);
328 G4LogicalVolume* l6a = new G4LogicalVolume(s6a, Materials::get("SUS304"), "l6a", 0, 0, 0);
329 l6a->SetVisAttributes(batt);
330 new G4PVPlacement(G4TranslateY3D(-45. / 2), l6a, "l6a_physical", l6, false, 0, overlap);
331
332 G4VSolid* s6b = new G4Box("leg_p6b", 60. / 2, 45. / 2, 160. / 2);
333 G4LogicalVolume* l6b = new G4LogicalVolume(s6b, Materials::get("SUS304"), "l6b", 0, 0, 0);
334 l6b->SetVisAttributes(batt);
335 double dy = 189. / 2 - 45 + 45. / 2;
336 new G4PVPlacement(G4TranslateY3D(dy), l6b, "l6b_physical", l6, false, 0, overlap);
337
338 G4VSolid* s6c = new G4Box("leg_p6c", 40. / 2, 45. / 2, 22.5 / 2);
339 G4LogicalVolume* l6c = new G4LogicalVolume(s6c, Materials::get("SUS304"), "l6c", 0, 0, 0);
340 l6c->SetVisAttributes(batt);
341 new G4PVPlacement(G4Translate3D(30 + 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 0, overlap);
342 new G4PVPlacement(G4Translate3D(30 + 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 1, overlap);
343 new G4PVPlacement(G4Translate3D(-30 - 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 2, overlap);
344 new G4PVPlacement(G4Translate3D(-30 - 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 3, overlap);
345
346 G4AssemblyVolume* support_leg = new G4AssemblyVolume();
347
348 support_leg->AddPlacedVolume(l1, t1);
349 // support_leg->AddPlacedVolume(l2,t2);
350 support_leg->AddPlacedVolume(l3, t3);
351 // support_leg->AddPlacedVolume(l4,t4);
352 // support_leg->AddPlacedVolume(l5,t5);
353 support_leg->AddPlacedVolume(l6, t6);
354
355 G4VSolid* s_all = new G4Box("leg_all", 140. / 2, 359. / 2, (257. - 5.) / 2);
356 G4LogicalVolume* l_all = new G4LogicalVolume(s_all, Materials::get("G4_AIR"), "l_all", 0, 0, 0);
357 l_all->SetVisAttributes(att("air"));
358 G4Transform3D tp = G4Translate3D(0, -359. / 2, -(257. - 5.) / 2);
359 support_leg->MakeImprint(l_all, tp, 0, overlap);
360
361
362 for (int i = 0; i < 8; i++)
363 new G4PVPlacement(gTrans * G4RotateX3D(M_PI)*G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0,
364 1496 - 185 + 359. / 2,
365 434 + 5 + (257. - 5.) / 2), l_all, "ECL_BackwardSupport_lall_physical", top, false, i, overlap);
366
367
368 G4VSolid* s1a = new G4Box("leg_p1a", 130. / 2, 178. / 2, 5. / 2);
369 G4LogicalVolume* l1a = new G4LogicalVolume(s1a, Materials::get("SUS304"), "l1a", 0, 0, 0);
370 l1a->SetVisAttributes(batt);
371 for (int i = 0; i < 8; i++)
372 new G4PVPlacement(gTrans * G4RotateX3D(M_PI)*G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0,
373 1496 - 185 + 178. / 2,
374 434 + 5 - 5. / 2), l1a, "l1a_physical", top, false, i, overlap);
375
376 delete support_leg;
377 }
378
379}
a Belle crystal in Geant4
Definition: BelleCrystal.h:37
BelleLathe class.
Definition: BelleLathe.h:67
double get_pa_box_height() const
Getter for preamplifier box height (hard-coded to be 2)
Definition: GeoECLCreator.h:88
const G4VisAttributes * att(const std::string &n) const
Define visual attributes.
G4LogicalVolume * get_preamp() const
Get Logical volume of preamplifier.
G4LogicalVolume * wrapped_crystal(const shape_t *s, const std::string &endcap, double wrapthickness)
Wrapped crystal.
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
double atan(double a)
atan for double
Definition: beamHelpers.h:34
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
const int c_NCrystalsForwardBarrel
Number of crystals in the forward and barrel ECL.
const std::vector< double > v3
MATLAB generated random vector.
struct for Point
Definition: BelleCrystal.h:31
placement struct
Definition: shapes.h:49
simple struct with z and r coordinates
Definition: BelleLathe.h:25

◆ barrel()

void barrel ( G4LogicalVolume &  _top)
private

Make the ECL barrel and then place elements inside it.

Definition at line 63 of file barrel.cc.

64{
65 G4LogicalVolume* top = &_top;
66
67 // vector<cplacement_t> bp = load_placements("/ecl/data/crystal_placement_barrel.dat");
68 vector<cplacement_t> bp = load_placements(m_sap, ECLParts::barrel);
69
70 vector<cplacement_t>::iterator fp = find_if(bp.begin(), bp.end(), [](const cplacement_t& p) {
71 const int ECL_barrel_part = 1001;
72 return p.nshape == ECL_barrel_part;
73 });
74 // global transformation before placing the whole forward part in the top logical volume
75 G4Transform3D gT = (fp == bp.end()) ? G4Translate3D(0, 0, 0) : get_transform(*fp);
76 if (fp != bp.end()) bp.erase(fp); // now not needed
77
78 const bool b_crystals = true;
79 const bool b_forward_support_legs = true;
80 const bool b_backward_support_legs = true;
81 const bool b_forward_support_ring = true;
82 const bool b_backward_support_ring = true;
83 const bool b_forward_support_wall = true;
84 const bool b_backward_support_wall = true;
85 const bool b_crystal_support = true;
86 const bool b_preamplifier = true;
87 const bool b_septum_walls = true;
88 const bool b_inner_support_wall = true;
89 const bool b_support_ribs = true;
90 const bool b_outer_plates = true;
91 const bool b_forward_part5 = true;
92 const bool b_crystal_holder = true;
93 const bool b_pipe = true;
94
95 const double phi0 = M_PI / 2 - M_PI / 72;
96 const double dphi = 2 * M_PI; //2*2*M_PI/72;
97
98 const int nseg = 72;
99
100 int overlap = m_overlap;
101
102 G4LogicalVolume* sectorlogical; // tilted sector
103 {
104 double R = 1251.6;
105 double zp = 4.9, zm = 2.5;
106 // zr_t bint[] = {{-1000.8-zm, 0}, {1982.6 + zp,0}, {2288.,198.2 - zp*tand(32.982)}, {2288.,404 + 12}, {-1219.,404 + 12}, {-1219., 404 -112.6-zm*tand(52.897)-3.2}};
107 zr_t bint[] = {{ -1000.8 - zm, 0}, {1982.6 + zp, 0}, {2288., 198.2 - zp * tand(32.982)}, {2288., 1582 - R}, {2380, 1582 - R}, {2380, 404 + 12}, { -1330., 404 + 12}, { -1330., 1543 + 84 - R}, { -1219., 1543 + 84 - R}, { -1219., 404 - 112.6 - zm * tand(52.897) - 3.2}};
108 for (unsigned int i = 0; i < sizeof(bint) / sizeof(zr_t); i++) bint[i].r += R;
109 std::vector<zr_t> contourb(bint, bint + sizeof(bint) / sizeof(zr_t));
110 G4VSolid* sect0 = new BelleLathe("sect0", phi0, 2 * M_PI, contourb);
111
112 G4Vector3D n0(cos(phi0), sin(phi0), 0);
113 // shift r0 by half of septum wall thickness so from one side of the sector the wall will be completely inside the sector volume
114 G4Point3D r0 = R * n0 + 0.5 / 2 * (G4RotateZ3D(M_PI / 2) * n0), rt0 = r0 + G4RotateZ3D(M_PI / 144 - M_PI / 72) * (G4Vector3D(0, 417,
115 0)), r1 = G4RotateZ3D(2 * M_PI / 72) * r0, rt1 = G4RotateZ3D(2 * M_PI / 72) * rt0;
116 G4Transform3D t0 = G4RotateZ3D(-M_PI / 144) * G4Translate3D(-r0);
117 G4Point3D tr0 = t0 * r0, trt0 = t0 * rt0, tr1 = t0 * r1, trt1 = t0 * rt1;
118 G4Vector3D tn1 = trt1 - tr1;
119 double v0 = -(tr1.y() - tr0.y()) / tn1.y();
120 double v1 = -(tr1.y() - trt0.y()) / tn1.y();
121
122 G4Point3D u0 = tr1 + v0 * tn1, u1 = tr1 + v1 * tn1;
123 G4VSolid* sect1 = new G4Trd("sect1", -u0.x() / 2, (trt0 - u1).x() / 2, 2000, 2000, trt0.y() / 2);
124 G4IntersectionSolid* sector = new G4IntersectionSolid("sector", sect0, sect1,
125 G4Translate3D(r0)*G4RotateZ3D(M_PI / 144)*G4Translate3D(u0.x() / 2, trt0.y() / 2, (2288 - 1219.) / 2)*G4RotateX3D(-M_PI / 2));
126 sectorlogical = new G4LogicalVolume(sector, Materials::get("G4_AIR"), "ECLBarrelSectorLogical", 0, 0, 0);
127 sectorlogical->SetVisAttributes(att("air"));
128
129 // Set up region for production cuts
130 G4Region* aRegion = new G4Region("ECLBarrelSector");
131 sectorlogical->SetRegion(aRegion);
132 aRegion->AddRootLogicalVolume(sectorlogical);
133
134 for (int i = 0; i < nseg; i++) {
135 double phi = i * M_PI / 36 - M_PI / 2;
136 new G4PVPlacement(gT * G4RotateZ3D(phi), sectorlogical, suf("ECLBarrelSectorPhysical", i), top, false, i, overlap);
137 }
138 }
139
140 // phi septum wall
141 if (b_septum_walls) {
142 double zh = (0.5 - 0.001) / 2;
143 G4ThreeVector psw[] = {G4ThreeVector(-1000.8, 0, -zh), G4ThreeVector(1982.6, 0, -zh), G4ThreeVector(2288, 198.2, -zh),
144 G4ThreeVector(2288, 404, -zh), G4ThreeVector(-1219., 404, -zh), G4ThreeVector(-1219., 404 - 112.6, -zh),
145 G4ThreeVector(-1000.8, 0, zh), G4ThreeVector(1982.6, 0, zh), G4ThreeVector(2288, 198.2, zh),
146 G4ThreeVector(2288, 404, zh), G4ThreeVector(-1219., 404, zh), G4ThreeVector(-1219., 404 - 112.6, zh)
147 };
148 G4VSolid* pswsolid = new BelleCrystal("pswsolid", 6, psw);
149 G4LogicalVolume* pswlogical = new G4LogicalVolume(pswsolid, Materials::get("A5052"), "pswlogical", 0, 0, 0);
150 pswlogical->SetVisAttributes(att("alum"));
151 new G4PVPlacement(G4RotateZ3D(2.5 / 180 * M_PI)*G4Translate3D(0, 1251.6 + 0.01,
152 0)*G4RotateZ3D(1.25 / 180 * M_PI)*G4RotateY3D(-M_PI / 2), pswlogical, "phiwallphysical", sectorlogical, false, 0, overlap);
153 }
154
155 if (b_forward_support_ring) {
156 zr_t vc4[] = {{0, 1449.2}, {0, 1582}, {90, 1582}, {90, 1582 - 75.7}};
157 std::vector<zr_t> contour4(vc4, vc4 + 4);
158 G4VSolid* part4solid = new BelleLathe("part4solid", phi0, dphi, contour4);
159 G4LogicalVolume* part4logical = new G4LogicalVolume(part4solid, Materials::get("SUS304"), "part4logical", 0, 0, 0);
160 part4logical->SetVisAttributes(att("iron"));
161 new G4PVPlacement(gT * G4Translate3D(0, 0, 2290), part4logical, "ECL_BarrelSupport_part4physical", top, false, 0, overlap);
162 }
163
164 if (b_backward_support_ring) {
165 zr_t vc6[] = {{0, 1543}, { -105, 1543}, { -105, 1543 + 84}, {0, 1543 + 84}};
166 std::vector<zr_t> contour6(vc6, vc6 + 4);
167 G4VSolid* part6solid = new BelleLathe("part6solid", phi0, dphi, contour6);
168 G4LogicalVolume* part6logical = new G4LogicalVolume(part6solid, Materials::get("SUS304"), "part6logical", 0, 0, 0);
169 part6logical->SetVisAttributes(att("iron"));
170 new G4PVPlacement(gT * G4Translate3D(0, 0, -1225), part6logical, "part6physical", top, false, 0, overlap);
171 }
172
173 if (b_forward_support_wall) {
174 zr_t vc2[] = {{0, 1180.1}, {0, 1249.5}, {54.3, 1249.5}, {445.8, 1503.6}, {445.8 + 4 * sind(32.982), 1503.6 - 4 * cosd(32.982)}, {54.3 + 4 * sind(32.982 / 2), 1249.5 - 4}, {4, 1249.5 - 4}, {4, 1180.1}};
175 std::vector<zr_t> contour2(vc2, vc2 + 8);
176 G4VSolid* part2solid = new BelleLathe("part2solid", phi0, dphi, contour2);
177 G4LogicalVolume* part2logical = new G4LogicalVolume(part2solid, Materials::get("A5083"), "part2logical", 0, 0, 0);
178 part2logical->SetVisAttributes(att("alum"));
179 new G4PVPlacement(gT * G4Translate3D(0, 0, 1930), part2logical, "ECL_BarrelSupport_part2physical", top, false, 0, overlap);
180 }
181
182 if (b_backward_support_wall) {
183 zr_t vc3[] = {{0, 1180.1}, { -4, 1180.1}, { -4, 1249.5 - 4}, { -61.9 - 4 * sind(52.897 / 2), 1249.5 - 4}, { -285, 1539.1 - 4 * sind(52.897 / 2)}, { -285, 1618.8}, { -285 + 4, 1618.8}, { -285 + 4, 1539.1}, { -61.9, 1249.5}, {0, 1249.5}};
184 std::vector<zr_t> contour3(vc3, vc3 + 10);
185 G4VSolid* part3solid = new BelleLathe("part3solid", phi0, dphi, contour3);
186 G4LogicalVolume* part3logical = new G4LogicalVolume(part3solid, Materials::get("A5083"), "part3logical", 0, 0, 0);
187 part3logical->SetVisAttributes(att("alum"));
188 new G4PVPlacement(gT * G4Translate3D(0, 0, -940), part3logical, "ECL_BarrelSupport_part3physical", top, false, 0, overlap);
189 }
190
191 if (b_inner_support_wall) {
192 G4VSolid* part1solid = new G4Tubs("part1solid", 1250.1, 1250.1 + 1.5, (988.5 + 1972.5) / 2, phi0, dphi);
193 G4LogicalVolume* part1logical = new G4LogicalVolume(part1solid, Materials::get("A5083"), "part1logical", 0, 0, 0);
194 part1logical->SetVisAttributes(att("alum2"));
195 new G4PVPlacement(gT * G4Translate3D(0, 0, (-988.5 + 1972.5) / 2), part1logical, "ECL_BarrelSupport_part1physical", top, false, 0,
196 overlap);
197 }
198
199 // cppcheck-suppress knownConditionTrueFalse
200 if (b_support_ribs || b_outer_plates || b_forward_part5) {
201 // double zh = (1219.+2288.)/2, dz = (-1219.+2288.)/2;
202 double zh = (1330. + 2380.) / 2, dz = (-1330. + 2380.) / 2;
203
204 if (b_support_ribs) {
205 G4ThreeVector vc9[] = {G4ThreeVector(0, 0, -zh), G4ThreeVector(25, 0, -zh), G4ThreeVector(25, 24.5 - 25.0 * tand(2.5), -zh), G4ThreeVector(0, 24.5, -zh),
206 G4ThreeVector(0, 0, zh), G4ThreeVector(25, 0, zh), G4ThreeVector(25, 24.5 - 25.0 * tand(2.5), zh), G4ThreeVector(0, 24.5, zh)
207 };
208 G4VSolid* part9solid = new BelleCrystal("part9solid", 4, vc9);
209 G4LogicalVolume* part9logical = new G4LogicalVolume(part9solid, Materials::get("SUS304"), "part9logical", 0, 0, 0);
210 part9logical->SetVisAttributes(att("silv"));
211 G4Transform3D tpart9a = G4RotateZ3D(2.5 / 180 * M_PI) * G4Translate3D(0, 1251.6,
212 0) * G4RotateZ3D(1.25 / 180 * M_PI) * G4Translate3D(0.5 / 2, 404 - 25 - 3.7 + 0.2, dz);
213 G4Transform3D tpart9b = G4RotateZ3D(-2.5 / 180 * M_PI) * G4Translate3D(0, 1251.6,
214 0) * G4RotateZ3D(1.25 / 180 * M_PI) * G4Translate3D(-0.5 / 2, 404 - 25 - 3.7 + 0.5, dz) * G4RotateY3D(M_PI);
215 new G4PVPlacement(tpart9a, part9logical, "part9aphysical", sectorlogical, false, 0, overlap);
216 new G4PVPlacement(tpart9b, part9logical, "part9bphysical", sectorlogical, false, 0, overlap);
217 }
218
219 if (b_outer_plates) {
220 G4VSolid* part10solid = new G4Box("part10solid", 143. / 2, 8. / 2, zh);
221 G4LogicalVolume* part10logical = new G4LogicalVolume(part10solid, Materials::get("SUS304"), "part10logical", 0, 0, 0);
222 part10logical->SetVisAttributes(att("iron2"));
223 G4Transform3D tpart10 = G4Translate3D(0, 1251.6, 0) * G4Translate3D(-404 * tand(1.25), 404 - 1.1,
224 dz) * G4RotateZ3D(0.3 / 180 * M_PI);
225 new G4PVPlacement(tpart10, part10logical, "part10physical", sectorlogical, false, 0, overlap);
226 }
227
228 // G4VSolid *part5solid = new G4Box("part5solid", 60./2, 45./2, 90./2);
229 // G4LogicalVolume* part5logical = new G4LogicalVolume(part5solid, Materials::get("SUS304"), "part5logical", 0, 0, 0);
230 // part5logical->SetVisAttributes(att("iron2"));
231 // new G4PVPlacement(G4RotateZ3D(2.796/180*M_PI)*G4Translate3D(0,1582+45./2,2380.-90./2), part5logical, "part5physical", sectorlogical, false, 1, overlap);
232
233
234 if (b_forward_part5) { // slice the element to fit into the sector
235 double R = 1251.6, phi1 = phi0 + 1.25 * M_PI / 180, phi2 = (2.796 + 90 - 5) / 180 * M_PI;
236 G4Vector3D n0(cos(phi0), sin(phi0), 0), n1(cos(phi1), sin(phi1), 0), n2(cos(phi2), sin(phi2), 0);
237 // shift r0 by half of septum wall thickness so from one side of the sector the wall will be completely inside the sector volume
238 G4Point3D r0 = R * n0 + 0.5 / 2 * (G4RotateZ3D(M_PI / 2) * n0);
239
240 // n2*(r0 + t*n1 - n2*1582) = 0;
241 double t0 = -n2 * (r0 - n2 * 1582) / (n2 * n1);
242 double t1 = -n2 * (r0 - n2 * (1582 + 45)) / (n2 * n1);
243
244 G4Point3D x0 = r0 + t0 * n1;
245 G4Point3D x1 = (t1 - t0) * n1;
246
247 Point_t c5a[] = {{0, 0}, {30 * cosd(5 - 2.796), -30 * sind(5 - 2.796)}, {30 * cosd(5 - 2.796) + 45 * cos(phi2), -30 * sind(5 - 2.796) + 45 * sin(phi2)}, {x1.x(), x1.y()}};
248 G4ThreeVector vc5a[8];
249 for (int i = 0; i < 4; i++) vc5a[i] = G4ThreeVector(c5a[i].x, c5a[i].y, -90. / 2);
250 for (int i = 0; i < 4; i++) vc5a[i + 4] = G4ThreeVector(c5a[i].x, c5a[i].y, 90. / 2);
251 G4VSolid* part5asolid = new BelleCrystal("part5asolid", 4, vc5a);
252 G4LogicalVolume* part5alogical = new G4LogicalVolume(part5asolid, Materials::get("SUS304"), "part5alogical", 0, 0, 0);
253 part5alogical->SetVisAttributes(att("iron2"));
254 new G4PVPlacement(G4RotateZ3D(2 * M_PI / 72)*G4Translate3D(x0.x(), x0.y(), 2380. - 90. / 2), part5alogical, "part5aphysical",
255 sectorlogical, false, 0, overlap);
256
257 Point_t c5b[] = {{ -30 * cosd(5 - 2.796), 30 * sind(5 - 2.796)}, {0, 0}, {x1.x(), x1.y()}, { -30 * cosd(5 - 2.796) + 44.8 * cos(phi2), 30 * sind(5 - 2.796) + 44.8 * sin(phi2)}};
258 G4ThreeVector vc5b[8];
259 for (int i = 0; i < 4; i++) vc5b[i] = G4ThreeVector(c5b[i].x, c5b[i].y, -90. / 2);
260 for (int i = 0; i < 4; i++) vc5b[i + 4] = G4ThreeVector(c5b[i].x, c5b[i].y, 90. / 2);
261 G4VSolid* part5bsolid = new BelleCrystal("part5bsolid", 4, vc5b);
262 G4LogicalVolume* part5blogical = new G4LogicalVolume(part5bsolid, Materials::get("SUS304"), "part5blogical", 0, 0, 0);
263 part5blogical->SetVisAttributes(att("iron2"));
264 new G4PVPlacement(G4Translate3D(x0.x(), x0.y(), 2380. - 90. / 2), part5blogical, "part5bphysical", sectorlogical, false, 0,
265 overlap);
266 }
267
268 }
269
270 if (b_forward_support_legs) {
271 G4ThreeVector vc11a[] = {G4ThreeVector(0, 0, -63. / 2), G4ThreeVector(35, 24.5, -63. / 2), G4ThreeVector(35, 158 - 35, -63. / 2), G4ThreeVector(0, 158 - 35, -63. / 2),
272 G4ThreeVector(0, 0, 63. / 2), G4ThreeVector(35, 24.5, 63. / 2), G4ThreeVector(35, 158 - 35, 63. / 2), G4ThreeVector(0, 158 - 35, 63. / 2)
273 };
274 G4VSolid* part11asolid = new BelleCrystal("part11asolid", 4, vc11a);
275 G4VSolid* part11bsolid = new G4Box("part11bsolid", 63. / 2, 35. / 2, 100. / 2);
276 G4LogicalVolume* part11alogical = new G4LogicalVolume(part11asolid, Materials::get("SUS304"), "part11alogical", 0, 0, 0);
277 part11alogical->SetVisAttributes(att("iron2"));
278 G4LogicalVolume* part11blogical = new G4LogicalVolume(part11bsolid, Materials::get("SUS304"), "part11blogical", 0, 0, 0);
279 part11blogical->SetVisAttributes(att("iron"));
280
281 double Ro = 158.4;
282 zr_t vsf[] = {{0, 0}, {35, 24.5}, {35, 158 - 35}, {100, 158 - 35}, {100, Ro}, {0, Ro}};
283 for (zr_t* it = vsf; it != vsf + 6; it++) {it->r += 1668 - 158;}
284 std::vector<zr_t> csf(vsf, vsf + 6);
285 G4VSolid* sf = new BelleLathe("sf", 0, 2 * M_PI, csf);
286 G4LogicalVolume* sfl = new G4LogicalVolume(sf, Materials::get("G4_AIR"), "supportfwd", 0, 0, 0);
287 sfl->SetVisAttributes(att("air"));
288 new G4PVPlacement(gT * G4Translate3D(0, 0, 2380), sfl, "supportfwdphysical", top, false, 0, overlap);
289
290 G4VSolid* sfs = new BelleLathe("sfs", -M_PI / 72, 2 * M_PI / 72, csf);
291 G4LogicalVolume* sfsl = new G4LogicalVolume(sfs, Materials::get("G4_AIR"), "supportfwdsector", 0, 0, 0);
292 sfsl->SetVisAttributes(att("air"));
293
294 new G4PVPlacement(G4Translate3D(1668 - 158, 0, 0)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2), part11alogical, "part11aphysical",
295 sfsl, false, 0, overlap);
296 new G4PVPlacement(G4Translate3D(1668 - 158, 0, 0)*G4RotateZ3D(-M_PI / 2)*G4Translate3D(0, (158. - 35.) + 35. / 2, 100. / 2),
297 part11blogical, "part11bphysical", sfsl, false, 0, overlap);
298 new G4PVReplica("supportfwdrep", sfsl, sfl, kPhi, nseg, 2 * M_PI / 72, (2.796 - 2.5)*M_PI / 180);
299 }
300
301 if (b_backward_support_legs) {
302 G4VSolid* part12asolid = new G4Box("part12asolid", 90. / 2, (118. - 35.) / 2, 35. / 2);
303 G4VSolid* part12bsolid = new G4Box("part12bsolid", 90. / 2, 35. / 2, 120. / 2);
304 G4LogicalVolume* part12alogical = new G4LogicalVolume(part12asolid, Materials::get("SUS304"), "part12alogical", 0, 0, 0);
305 part12alogical->SetVisAttributes(att("iron2"));
306 G4LogicalVolume* part12blogical = new G4LogicalVolume(part12bsolid, Materials::get("SUS304"), "part12blogical", 0, 0, 0);
307 part12blogical->SetVisAttributes(att("iron"));
308
309 double Rb = 118.7;
310 zr_t vsb[] = {{0, 0}, {0, Rb}, { -120, Rb}, { -120, 118 - 35}, { -35, 118 - 35}, { -35, 0}};
311 for (zr_t* it = vsb; it != vsb + 6; it++) {it->r += 1668 - 118;}
312 std::vector<zr_t> csb(vsb, vsb + 6);
313 G4VSolid* sb = new BelleLathe("sb", 0, 2 * M_PI, csb);
314 G4LogicalVolume* sbl = new G4LogicalVolume(sb, Materials::get("G4_AIR"), "supportbkw", 0, 0, 0);
315 sbl->SetVisAttributes(att("air"));
316 new G4PVPlacement(gT * G4Translate3D(0, 0, -1330), sbl, "supportbkwphysical", top, false, 0, overlap);
317
318 G4VSolid* sbs = new BelleLathe("sbs", -M_PI / 72, 2 * M_PI / 72, csb);
319 G4LogicalVolume* sbsl = new G4LogicalVolume(sbs, Materials::get("G4_AIR"), "supportbkwsector", 0, 0, 0);
320 sbsl->SetVisAttributes(att("air"));
321
322 new G4PVPlacement(G4Translate3D(1668 - 118, 0, 0)*G4RotateZ3D(-M_PI / 2)*G4Translate3D(0, (118. - 35.) / 2, -35. / 2),
323 part12alogical, "part12aphysical", sbsl, false, 0, overlap);
324 new G4PVPlacement(G4Translate3D(1668 - 118, 0, 0)*G4RotateZ3D(-M_PI / 2)*G4Translate3D(0, (118. - 35.) + 35. / 2, -120. / 2),
325 part12blogical, "part12bphysical", sbsl, false, 0, overlap);
326
327 new G4PVReplica("supportbkwrep", sbsl, sbl, kPhi, nseg, 2 * M_PI / 72, (2.796 - 2.5)*M_PI / 180);
328 }
329
330 if (b_septum_walls) {
331 struct zwall_t { double zlow, zup;};
332 zwall_t zs[] = {{ -682.1, -883.9}, { -445.9, -571.7}, { -220.8, -274.2}, {0, 0}, {220.8, 274.2}, {445.9, 571.7}, {682.1, 883.9}, {936.6, 1220.5}, {1217.1, 1591.2}, {1532., 2007.4}};
333
334 const double dz = 0.5 / 2; // half of theta fin thickness
335
336 G4Vector3D n0(sind(90 + 2.5 + 1.25), -cosd(90 + 2.5 + 1.25), 0); // normal vector of the left sector wall pointing inside the sector
337 G4Point3D r0 = G4Point3D(-1251.6 * sind(2.5), 0, 0) + dz * n0; // point in the left sector wall shifted by phi fin
338 G4Vector3D n1(-sind(90 - 2.5 + 1.25), cosd(90 - 2.5 + 1.25),
339 0); // normal vector of the right sector wall pointing inside the sector
340 G4Point3D r1 = G4Point3D(1251.6 * sind(2.5), 0, 0) + dz * n1; // point in the right sector wall shifted by phi fin
341
342 G4Vector3D n2(0, 1, 0); G4Point3D r2(0, 0, 0), r3(0, 310, 0);
343
344 for (unsigned int i = 0; i < sizeof(zs) / sizeof(zwall_t); i++) {
345 double th = atan2(404, -(zs[i].zlow - zs[i].zup));
346
347 G4Vector3D nf(0, -cos(th), sin(th));
348 // take into account finite thickness of the theta fin
349 G4Point3D rf(0, 0, copysign(dz / sin(th), cos(th)));
350
351 // calculate three plane intersection point
352 auto inter = [](const G4Vector3D & l_n0, const G4Point3D & l_r0,
353 const G4Vector3D & l_n1, const G4Point3D & l_r1,
354 const G4Vector3D & l_n2, const G4Point3D & l_r2) -> G4Point3D {
355 CLHEP::HepMatrix A(3, 3);
356 CLHEP::HepVector B(3);
357 A[0][0] = l_n0.x(), A[0][1] = l_n0.y(), A[0][2] = l_n0.z();
358 A[1][0] = l_n1.x(), A[1][1] = l_n1.y(), A[1][2] = l_n1.z();
359 A[2][0] = l_n2.x(), A[2][1] = l_n2.y(), A[2][2] = l_n2.z();
360
361 B[0] = l_r0 * l_n0;
362 B[1] = l_r1 * l_n1;
363 B[2] = l_r2 * l_n2;
364
365 CLHEP::HepVector r = A.inverse() * B;
366 G4Point3D res(r[0], r[1], r[2]);
367 return res;
368 };
369
370 G4Transform3D tfin(G4RotateX3D(-th + M_PI / 2)), tfini(G4RotateX3D(th - M_PI / 2));
371 G4ThreeVector t0 = tfini * inter(n0, r0, n2, r2, nf, rf);
372 G4ThreeVector t1 = tfini * inter(n1, r1, n2, r2, nf, rf);
373 G4ThreeVector t2 = tfini * inter(n0, r0, n2, r3, nf, rf);
374 G4ThreeVector t3 = tfini * inter(n1, r1, n2, r3, nf, rf);
375 // G4cout<<t0<<" "<<t1<<" "<<t2<<" "<<t3<<" "<<G4RotateX3D(th-M_PI/2)*G4Point3D(0,404,zs[i].zup-zs[i].zlow)<<G4endl;
376
377 G4ThreeVector thfin[8];
378 thfin[0] = G4ThreeVector(t0.x(), t0.y(), -dz);
379 thfin[1] = G4ThreeVector(t1.x(), t1.y(), -dz);
380 thfin[2] = G4ThreeVector(t3.x(), t3.y(), -dz);
381 thfin[3] = G4ThreeVector(t2.x(), t2.y(), -dz);
382 thfin[4] = G4ThreeVector(t0.x(), t0.y(), +dz);
383 thfin[5] = G4ThreeVector(t1.x(), t1.y(), +dz);
384 thfin[6] = G4ThreeVector(t3.x(), t3.y(), +dz);
385 thfin[7] = G4ThreeVector(t2.x(), t2.y(), +dz);
386
387 G4VSolid* thfinsolid = new BelleCrystal("thfinsolid", 4, thfin);
388 G4LogicalVolume* thfinlogical = new G4LogicalVolume(thfinsolid, Materials::get("A5052"), "thfinlogical", 0, 0, 0);
389 thfinlogical->SetVisAttributes(att("iron2"));
390 new G4PVPlacement(G4Translate3D(0, 1251.6, zs[i].zlow)*tfin, thfinlogical, "thfinphysical", sectorlogical, false, 0, overlap);
391 }
392 }
393
394 if (b_pipe) {
395 double zh = (1219. + 2288.) / 2, dz = (-1219. + 2288.) / 2;
396 G4VSolid* tubec = new G4Tubs("tubec", 0, 12. / 2, zh, 0, 2 * M_PI);
397 G4VSolid* tubew = new G4Tubs("tubew", 0, 10. / 2, zh, 0, 2 * M_PI);
398 G4LogicalVolume* tubeclogical = new G4LogicalVolume(tubec, Materials::get("C1220"), "tubec", 0, 0, 0);
399 G4LogicalVolume* tubewlogical = new G4LogicalVolume(tubew, Materials::get("G4_WATER"), "tubew", 0, 0, 0);
400 new G4PVPlacement(G4Translate3D(0, 0, 0), tubewlogical, "tubewphysical", tubeclogical, false, 1, overlap);
401 tubeclogical->SetVisAttributes(att("iron2"));
402 new G4PVPlacement(G4Translate3D(25, 1640, dz), tubeclogical, "tubecphysical", sectorlogical, false, 0, overlap);
403 new G4PVPlacement(G4Translate3D(-35, 1640, dz), tubeclogical, "tubecphysical", sectorlogical, false, 1, overlap);
404 }
405
406 if (b_preamplifier) {
407 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
408 G4Transform3D twc = get_transform(*it);
409 int indx = it - bp.begin();
410 auto pv = new G4PVPlacement(twc * G4TranslateZ3D(300 / 2 + 0.250 + get_pa_box_height() / 2), get_preamp(),
411 suf("phys_barrel_preamplifier", indx),
412 sectorlogical, false, 72 + 9 * (indx / 2) + (indx % 2), 0);
413 if (overlap)pv->CheckOverlaps(1000);
414 }
415 }
416
417 if (b_crystal_holder) {
418 G4VSolid* holder0 = new G4Box("holder0", 19. / 2, 38. / 2, 6. / 2);
419 G4VSolid* holder1 = new G4Box("holder1", 7. / 2, 26. / 2, 7. / 2);
420 G4VSolid* holder = new G4SubtractionSolid("holder", holder0, holder1);
421 G4LogicalVolume* holderlogical = new G4LogicalVolume(holder, Materials::get("A5052"), "holderlogical", 0, 0, 0);
422 holderlogical->SetVisAttributes(att("holder"));
423 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
424 G4Transform3D twc = get_transform(*it);
425 int indx = it - bp.begin();
426 auto pv = new G4PVPlacement(twc * G4TranslateZ3D(300 / 2 + 0.250 + get_pa_box_height() + 38. / 2)*G4RotateZ3D(M_PI / 2)*G4RotateX3D(
427 M_PI / 2), holderlogical, "holderphysical", sectorlogical, false, indx, 0);
428 if (overlap)pv->CheckOverlaps(1000);
429 }
430 }
431
432 if (b_crystals) {
433 // vector<shape_t*> cryst = load_shapes("/ecl/data/crystal_shape_barrel.dat");
434 vector<shape_t*> cryst = load_shapes(m_sap, ECLParts::barrel);
435 vector<G4LogicalVolume*> wrapped_crystals;
436 for (auto it = cryst.begin(); it != cryst.end(); it++) {
437 shape_t* s = *it;
438 wrapped_crystals.push_back(wrapped_crystal(s, "barrel", 0.17 - 0.0325));
439 }
440
441 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
442 const cplacement_t& t = *it;
443 auto s = find_if(cryst.begin(), cryst.end(), [&t](const shape_t* shape) {return shape->nshape == t.nshape;});
444 if (s == cryst.end()) continue;
445
446 G4Transform3D twc = get_transform(t);
447 int indx = it - bp.begin();
448 new G4PVPlacement(twc, wrapped_crystals[s - cryst.begin()], suf("ECLBarrelWrappedCrystal_Physical", indx), sectorlogical, false,
449 72 + 9 * (indx / 2) + (indx % 2), overlap);
450 }
451
452 for (shape_t* shape : cryst) {
453 delete shape;
454 }
455 }
456
457 if (b_crystal_support) {
458 G4LogicalVolume* tl[8];
459 for (int i = 0; i < 8; i++) {
460 G4VSolid* t = get_crystal_support(i);
461 tl[i] = new G4LogicalVolume(t, Materials::get("SUS304"), "crystal_support_logical", 0, 0, 0);
462 tl[i]->SetVisAttributes(att("iron2"));
463 }
464
465 double offset[] = { -945, -607, -305, -45, 45, 305, 607, 909, 1240, 1628, 2040};
466 const int ns[] = {8, 5, 6, 7, 7, 6, 5, 4, 3, 2, 1};
467 const double offsetb[] = {263.4 - 210 - 30, 27.4, 10.6, 16.1, 26.2, 21.7, 35, 22 + 240};
468 bool flip[] = {1, 0, 0, 0, 1, 1, 1, 1, 1, 1, 1};
469 for (int i = 0; i < 11; i++) {
470 G4Transform3D tflip = flip[i] ? G4RotateY3D(-M_PI / 2) : G4RotateY3D(M_PI / 2);
471 G4Transform3D tr0(G4RotateZ3D(2.5 / 180 * M_PI)*G4Translate3D(0, 1251.6,
472 0)*G4RotateZ3D(1.25 / 180 * M_PI)*G4Translate3D(25 + 3. / 2 + 0.5 / 2, 404 - 25 - 3.7 + 24.5 - 25.0 * tand(2.5) - 0.5,
473 offset[i])*tflip * G4Translate3D(-offsetb[ns[i] - 1], 0, 0)*G4RotateX3D(M_PI));
474 new G4PVPlacement(tr0, tl[ns[i] - 1], "crystal_support_physical", sectorlogical, false, 0, overlap);
475 G4Transform3D tr1(G4RotateZ3D(-2.5 / 180 * M_PI)*G4Translate3D(0, 1251.6,
476 0)*G4RotateZ3D(1.25 / 180 * M_PI)*G4Translate3D(-25 - 3. / 2 - 0.5 / 2, 404 - 25 - 3.7 + 24.5 - 25.0 * tand(2.5) - 0.5,
477 offset[i])*tflip * G4Translate3D(-offsetb[ns[i] - 1], 0, 0)*G4RotateX3D(M_PI));
478 new G4PVPlacement(tr1, tl[ns[i] - 1], "crystal_support_physical", sectorlogical, false, 1, overlap);
479 }
480 }
481
482}
double R
typedef autogenerated by FFTW
const std::vector< double > v1
MATLAB generated random vector.
double r
r coordinate
Definition: BelleLathe.h:27

◆ create()

void create ( const GearDir content,
G4LogicalVolume &  topVolume,
geometry::GeometryTypes  type 
)
overridevirtual

Function to actually create the geometry, has to be overridden by derived classes.

Parameters
contentGearDir pointing to the parameters which should be used for construction
topVolumeTop volume in which the geometry has to be placed
typeType of geometry to be build

Implements CreatorBase.

Definition at line 77 of file GeoECLCreator.cc.

78{
79 m_sensediode = new SensitiveDiode("ECLSensitiveDiode");
80 G4SDManager::GetSDMpointer()->AddNewDetector(m_sensediode);
81
82 ECLCrystalsShapeAndPosition crystals = loadCrystalsShapeAndPosition();
83 m_sap = &crystals;
84
85 forward(topVolume);
86 barrel(topVolume);
87 backward(topVolume);
88}
void barrel(G4LogicalVolume &)
Make the ECL barrel and then place elements inside it.
Definition: barrel.cc:63
void backward(G4LogicalVolume &)
Place elements inside the backward endcap.
Definition: backward.cc:42
void forward(G4LogicalVolume &)
Place elements inside the forward endcap.
Definition: forward.cc:44
Class for ECL Sensitive Detector.

◆ createFromDB()

void createFromDB ( const std::string &  name,
G4LogicalVolume &  topVolume,
geometry::GeometryTypes  type 
)
overridevirtual

Function to create the geometry from the Database.

Parameters
namename of the component in the database, could be used to disambiguate multiple components created with the same creator
topVolumeTop volume in which the geometry has to be placed
typeType of geometry to be build

Reimplemented from CreatorBase.

Definition at line 57 of file GeoECLCreator.cc.

58{
59 m_sensediode = new SensitiveDiode("ECLSensitiveDiode");
60 G4SDManager::GetSDMpointer()->AddNewDetector(m_sensediode);
61
63 if (!crystals.isValid()) B2FATAL("No crystal's data in the database.");
64 m_sap = &(*crystals);
65
66 forward(topVolume);
67 barrel(topVolume);
68 backward(topVolume);
69}
bool isValid() const
Check whether a valid object was obtained from the database.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21

◆ createPayloads()

void createPayloads ( const GearDir content,
const IntervalOfValidity iov 
)
overridevirtual

Function to create the geometry database.

This function should be implemented to convert Gearbox parameters to one ore more database payloads

Parameters
contentGearDir pointing to the parameters which should be used for construction
iovinterval of validity to use when generating payloads

Reimplemented from CreatorBase.

Definition at line 71 of file GeoECLCreator.cc.

72{
73 ECLCrystalsShapeAndPosition crystals = loadCrystalsShapeAndPosition();
75}
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42
bool storeData(const std::string &name, TObject *object, const IntervalOfValidity &iov)
Store an object in the database.
Definition: Database.cc:141

◆ defineVisAttributes()

void defineVisAttributes ( )
private

Define visual attributes.

Definition at line 119 of file GeoECLCreator.cc.

120{
121 m_atts["wrap"] = new G4VisAttributes(G4Colour(0.5, 0.5, 1.0));
122 m_atts["cryst"] = new G4VisAttributes(G4Colour(0.7, 0.7, 1.0));
123
124 m_atts["iron"] = new G4VisAttributes(G4Colour(1., 0.1, 0.1));
125 m_atts["iron2"] = new G4VisAttributes(G4Colour(1., 0.5, 0.5));
126 m_atts["alum"] = new G4VisAttributes(G4Colour(0.25, 0.25, 1.0, 0.5));
127 m_atts["alum2"] = new G4VisAttributes(G4Colour(0.5, 0.5, 1.0));
128 m_atts["silv"] = new G4VisAttributes(G4Colour(0.9, 0., 0.9));
129 m_atts["air"] = new G4VisAttributes(G4Colour(1., 1., 1.));
130 m_atts["air"]->SetVisibility(false);
131 m_atts["preamp"] = new G4VisAttributes(G4Colour(0.1, 0.1, 0.8));
132 m_atts["plate"] = new G4VisAttributes(G4Colour(0.2, 0.8, 0.2));
133 m_atts["connector"] = new G4VisAttributes(G4Colour(0.1, 0.1, 0.1));
134 m_atts["capacitor"] = new G4VisAttributes(G4Colour(0.1, 0.1, 0.8));
135 m_atts["holder"] = new G4VisAttributes(G4Colour(0.4, 0.8, 0.8));
136 m_atts["asolid6"] = new G4VisAttributes(G4Colour(1., 0.3, 0.2));
137 m_atts["asolid7"] = new G4VisAttributes(G4Colour(1., 0.3, 0.2));
138 m_atts["asolid13"] = new G4VisAttributes(G4Colour(1., 0.5, 0.5));
139}

◆ forward()

void forward ( G4LogicalVolume &  _top)
private

Place elements inside the forward endcap.

Definition at line 44 of file forward.cc.

45{
46 G4LogicalVolume* top = &_top;
47
48 const bool sec = 0;
49 const double phi0 = 0;
50 const double dphi = sec ? M_PI / 16 : 2 * M_PI;
51
52 const bool b_inner_support_ring = true;
53 const bool b_outer_support_ring = true;
54 const bool b_support_wall = true;
55 const bool b_ribs = true;
56 const bool b_septum_wall = true;
57 const bool b_crystals = true;
58 const bool b_preamplifier = true;
59 const bool b_support_leg = true;
60 const bool b_support_structure_13 = true;
61 const bool b_support_structure_15 = true;
62 bool b_connectors = true;
63 bool b_boards = true;
64 const bool b_cover = true;
65
66 b_connectors &= b_support_structure_15;
67 b_boards &= b_support_structure_15;
68
69 int overlap = m_overlap;
70
71 // vector<cplacement_t> bp = load_placements("/ecl/data/crystal_placement_forward.dat");
72 vector<cplacement_t> bp = load_placements(m_sap, ECLParts::forward);
73 vector<cplacement_t>::iterator fp = find_if(bp.begin(), bp.end(), [](const cplacement_t& p) {
74 const int ECL_forward_part = 1000;
75 return p.nshape == ECL_forward_part;
76 });
77 // global transformation before placing the whole forward part in the top logical volume
78 G4Transform3D gTrans = (fp == bp.end()) ? G4Translate3D(0, 0, 1960) : get_transform(*fp);
79 // since the calorimeter supporting structures attach to the yoke at
80 // the same place we need to modify supporting legs by milling
81 // necessary thickness of steel when moving forward part outward in
82 // z-direction
83 double milled_thickness = gTrans.dz() - 1960; // [mm]
84 if (fp != bp.end()) bp.erase(fp); // now not needed
85
86 const double th0 = 13.12, th1 = 32.98;
87 const double ZT = 437, ZI = 434, RI = 431, RIp = 532.2, RC = 1200.4, RT = 1415;
88 const double thinnerPart_translation = 1.95; // translation of the thin inner ring of the support ring
89
90 if (b_inner_support_ring) {
91 // This object has a (newly) translated part to eliminate overlaps with ARICH - 2021-07.
92 zr_t vc1[] = {{ZI - 487, 410 - thinnerPart_translation}, {ZT - (RIp - 410 - 20 / cosd(th0)) / tand(th0), 410 - thinnerPart_translation}, {ZT - (RIp - 410 - 20 / cosd(th0)) / tand(th0), 410}, {ZT, RIp - 20 / cosd(th0)}, {ZT, RIp}, {3., RI}, {3., 418 - thinnerPart_translation}, {ZI - 487, 418 - thinnerPart_translation}};
93 std::vector<zr_t> contour1(vc1, vc1 + sizeof(vc1) / sizeof(zr_t));
94 G4VSolid* part1solid = new BelleLathe("fwd_part1solid", phi0, dphi, contour1);
95 G4LogicalVolume* part1logical = new G4LogicalVolume(part1solid, Materials::get("SUS304"), "part1logical", 0, 0, 0);
96 part1logical->SetVisAttributes(att("iron"));
97 new G4PVPlacement(gTrans, part1logical, "ECL_ForwardSupport_part1physical", top, false, 0, overlap);
98 }
99
100 if (b_support_wall) {
101 // solving equation to get L : 3+L*cosd(th1)+1.6*cosd(th1+90) = 434 + 3 - 107.24
102 double L = (ZT - 107.24 - 3 - 1.6 * cosd(th1 + 90)) / cosd(th1);
103 double R0 = 418, R1 = RC;
104 zr_t vc23[] = {{0, R0}, {3, R0}, {3, R1}, {3 + L * cosd(th1), R1 + L * sind(th1)},
105 {3 + L * cosd(th1) + 1.6 * cosd(th1 + 90), R1 + L * sind(th1) + 1.6 * sind(th1 + 90)}, {3 + 1.6 * cosd(th1 + 90), R1 + 1.6 * sind(th1 + 90)}, {0, R1}
106 };
107 std::vector<zr_t> contour23(vc23, vc23 + sizeof(vc23) / sizeof(zr_t));
108 G4VSolid* part23solid = new BelleLathe("fwd_part23solid", phi0, dphi, contour23);
109 G4LogicalVolume* part23logical = new G4LogicalVolume(part23solid, Materials::get("A5052"), "part23logical", 0, 0, 0);
110 part23logical->SetVisAttributes(att("alum"));
111 new G4PVPlacement(gTrans, part23logical, "ECL_ForwardSupport_part23physical", top, false, 0, overlap);
112 }
113
114 if (b_outer_support_ring) {
115 zr_t vc4[] = {{3 + (RT - 20 - RC) / tand(th1), RT - 20}, {ZT, RT - 20}, {ZT, RT}, {3 + (RT - RC) / tand(th1), RT}};
116 std::vector<zr_t> contour4(vc4, vc4 + sizeof(vc4) / sizeof(zr_t));
117 G4VSolid* part4solid = new BelleLathe("fwd_part4solid", phi0, dphi, contour4);
118 G4LogicalVolume* part4logical = new G4LogicalVolume(part4solid, Materials::get("SUS304"), "part4logical", 0, 0, 0);
119 part4logical->SetVisAttributes(att("iron"));
120 new G4PVPlacement(gTrans, part4logical, "ECL_ForwardSupport_part4physical", top, false, 0, overlap);
121 }
122
123 zr_t cont_array_in[] = {{3., RI}, {ZT, RIp}, {ZT, RT - 20}, {3 + (RT - 20 - RC) / tand(th1), RT - 20}, {3, RC}};
124 std::vector<zr_t> contour_in(cont_array_in, cont_array_in + sizeof(cont_array_in) / sizeof(zr_t));
125 G4VSolid* innervolume_solid = new BelleLathe("fwd_innervolume_solid", 0, 2 * M_PI, contour_in);
126 G4LogicalVolume* innervolume_logical = new G4LogicalVolume(innervolume_solid, Materials::get("G4_AIR"),
127 "innervolume_logical", 0, 0, 0);
128 innervolume_logical->SetVisAttributes(att("air"));
129
130 // Set up region for production cuts
131 G4Region* aRegion = new G4Region("ECLForwardEnvelope");
132 innervolume_logical->SetRegion(aRegion);
133 aRegion->AddRootLogicalVolume(innervolume_logical);
134
135 new G4PVPlacement(gTrans, innervolume_logical, "ECLForwardPhysical", top, false, 0, overlap);
136
137 G4VSolid* innervolumesector_solid = new BelleLathe("fwd_innervolumesector_solid", -M_PI / 8, M_PI / 4, contour_in);
138 G4LogicalVolume* innervolumesector_logical = new G4LogicalVolume(innervolumesector_solid, Materials::get("G4_AIR"),
139 "innervolumesector_logical", 0, 0, 0);
140 innervolumesector_logical->SetVisAttributes(att("air"));
141 new G4PVReplica("ECLForwardSectorPhysical", innervolumesector_logical, innervolume_logical, kPhi, 8, M_PI / 4, 0);
142
143 if (b_ribs) {
144 double H = 60, W = 20;
145 double X0 = RIp, X1 = RT - 20;
146 G4TwoVector r0o(X1, 0), r1o(X1 * sqrt(1 - pow(W / X1, 2)), W);
147 double beta = asin(W / X0);
148 G4TwoVector r0i(X0 / cos(beta / 2), 0), r1i(X0 * cos(beta / 2) - tan(beta / 2) * (W - X0 * sin(beta / 2)), W);
149 double dxymzp = (r0o - r0i).x(), dxypzp = (r1o - r1i).x();
150 double theta = atan(tand(th0) / 2);
151 double dxymzm = dxymzp + tand(th0) * H, dxypzm = dxypzp + tand(th0) * H;
152
153 G4TwoVector m0 = (r0i + r0o) * 0.5, m1 = (r1i + r1o) * 0.5, dm = m1 - m0;
154 double alpha = atan(dm.x() / dm.y());
155
156 G4VSolid* solid6_p1 = new G4Trap("fwd_solid6_p1", H / 2, theta, 0, W / 2, dxymzm / 2, dxypzm / 2, alpha, W / 2, dxymzp / 2,
157 dxypzp / 2,
158 alpha);
159 G4LogicalVolume* lsolid6_p1 = new G4LogicalVolume(solid6_p1, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
160 lsolid6_p1->SetVisAttributes(att("iron"));
161 G4Transform3D tsolid6_p1(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, W / 2, ZT - H / 2));
162 new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*tsolid6_p1, lsolid6_p1, "psolid6_p1", innervolumesector_logical, false, 0, overlap);
163 new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p1, lsolid6_p1, "psolid6_p2", innervolumesector_logical, false, 0, overlap);
164
165 H = 40;
166 dxymzm = dxymzp + tand(th0) * H, dxypzm = dxypzp + tand(th0) * H;
167 G4VSolid* solid6_p2 = new G4Trap("fwd_solid6_p2", H / 2, theta, 0, W / 2, dxypzm / 2, dxymzm / 2, -alpha, W / 2, dxypzp / 2,
168 dxymzp / 2,
169 -alpha);
170 G4LogicalVolume* lsolid6_p2 = new G4LogicalVolume(solid6_p2, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
171 lsolid6_p2->SetVisAttributes(att("iron"));
172 G4Transform3D tsolid6_p2(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, -W / 2, ZT - H / 2));
173 new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p2, lsolid6_p2, "psolid6_p3", innervolumesector_logical, false, 0, overlap);
174 new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid6_p2, lsolid6_p2, "psolid6_p4", innervolumesector_logical, false, 0, overlap);
175
176 double hpad = 148.4;
177 G4VSolid* solid7_p8 = new G4Box("fwd_solid7_p8", hpad / 2, (140. - 40) / 2 / 2, 40. / 2);
178 G4LogicalVolume* lsolid7 = new G4LogicalVolume(solid7_p8, Materials::get("SUS304"), "lsolid7", 0, 0, 0);
179 lsolid7->SetVisAttributes(att("iron"));
180 double dx = sqrt(X1 * X1 - 70 * 70) - hpad / 2;
181 G4Transform3D tsolid7_p1(G4Translate3D(dx, -20 - 25, ZT - 40. / 2));
182 new G4PVPlacement(tsolid7_p1, lsolid7, "psolid7_p1", innervolumesector_logical, false, 0, overlap);
183 G4Transform3D tsolid7_p2(G4Translate3D(dx, 20 + 25, ZT - 40. / 2));
184 new G4PVPlacement(tsolid7_p2, lsolid7, "psolid7_p2", innervolumesector_logical, false, 0, overlap);
185
186 double L = X1 - (X0 - tand(th0) * 40) - 10;
187 G4VSolid* solid13 = new G4Box("fwd_solid13", L / 2, 5. / 2, 18. / 2);
188 G4LogicalVolume* lsolid13 = new G4LogicalVolume(solid13, Materials::get("SUS304"), "lsolid13", 0, 0, 0);
189 lsolid13->SetVisAttributes(att("iron"));
190 G4Transform3D tsolid13(G4TranslateZ3D(ZT - 60 + 18. / 2)*G4TranslateY3D(-5. / 2 - 0.5 / 2)*G4TranslateX3D(X0 - tand(
191 th0) * 40 + L / 2 + 5));
192 new G4PVPlacement(tsolid13, lsolid13, "psolid13_p1", innervolumesector_logical, false, 0, overlap);
193 new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid13, lsolid13, "psolid13_p2", innervolumesector_logical, false, 0, overlap);
194 }
195
196
197 double zsep = 125;
198 if (b_septum_wall) {
199 double d = 5;
200 Point_t vin[] = {{ZT - zsep, RIp - tand(th0)* zsep}, {ZT - 60, RIp - tand(th0) * 60}, {ZT - 60, RT - 20 - d}, {ZT - zsep, RT - 20 - d}};
201 const int n = sizeof(vin) / sizeof(Point_t);
202 Point_t c = centerofgravity(vin, vin + n);
203 G4ThreeVector contour_swall[n * 2];
204 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2);
205 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2);
206
207 G4VSolid* septumwall_solid = new BelleCrystal("fwd_septumwall_solid", n, contour_swall);
208
209 G4LogicalVolume* septumwall_logical = new G4LogicalVolume(septumwall_solid, Materials::get("A5052"),
210 "septumwall_logical", 0, 0, 0);
211 septumwall_logical->SetVisAttributes(att("alum2"));
212 new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0), septumwall_logical,
213 "septumwall_physical", innervolumesector_logical, false, 0, overlap);
214
215 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
216 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
217
218 G4VSolid* septumwall2_solid = new BelleCrystal("fwd_septumwall2_solid", n, contour_swall);
219
220 G4LogicalVolume* septumwall2_logical = new G4LogicalVolume(septumwall2_solid, Materials::get("A5052"),
221 "septumwall2_logical", 0, 0, 0);
222 septumwall2_logical->SetVisAttributes(att("alum2"));
223 new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0.5 / 2 / 2),
224 septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 0, overlap);
225 new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, -0.5 / 2 / 2),
226 septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 1, overlap);
227 }
228
229 zr_t vcr[] = {{3., RI}, {ZT - zsep, RIp - tand(th0)* zsep}, {ZT - zsep, RT - 20}, {3 + (RT - 20 - RC) / tand(th1), RT - 20}, {3, RC}};
230 std::vector<zr_t> ccr(vcr, vcr + sizeof(vcr) / sizeof(zr_t));
231 G4VSolid* crystalvolume_solid = new BelleLathe("fwd_crystalvolume_solid", 0, M_PI / 8, ccr);
232 G4LogicalVolume* crystalvolume_logical = new G4LogicalVolume(crystalvolume_solid, Materials::get("G4_AIR"),
233 "crystalvolume_logical", 0, 0, 0);
234 crystalvolume_logical->SetVisAttributes(att("air"));
235 new G4PVPlacement(G4RotateZ3D(-M_PI / 8), crystalvolume_logical, "ECLForwardCrystalSectorPhysical_0", innervolumesector_logical,
236 false, 0, overlap);
237 new G4PVPlacement(G4RotateZ3D(0), crystalvolume_logical, "ECLForwardCrystalSectorPhysical_1", innervolumesector_logical, false, 1,
238 overlap);
239
240 if (b_septum_wall) {
241 double d = 5, aRC = RC - 30e-6;
242 Point_t vin[] = {{3., RI}, {ZT - zsep, RIp - tand(th0)* zsep}, {ZT - zsep, RT - 20 - d}, {3 + (RT - 20 - d - aRC) / tand(th1), RT - 20 - d}, {3, aRC}};
243 const int n = sizeof(vin) / sizeof(Point_t);
244 Point_t c = centerofgravity(vin, vin + n);
245 G4ThreeVector contour_swall[n * 2];
246
247 for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
248 for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
249
250 G4VSolid* septumwall3_solid = new BelleCrystal("fwd_septumwall3_solid", n, contour_swall);
251
252 G4LogicalVolume* septumwall3_logical = new G4LogicalVolume(septumwall3_solid, Materials::get("A5052"),
253 "septumwall3_logical", 0, 0, 0);
254 septumwall3_logical->SetVisAttributes(att("alum2"));
255 new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0.5 / 2 / 2), septumwall3_logical,
256 "septumwall3_physical_0", crystalvolume_logical, false, 0, overlap);
257 new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, -0.5 / 2 / 2),
258 septumwall3_logical, "septumwall3_physical_1", crystalvolume_logical, false, 1, overlap);
259 }
260
261 if (b_crystals) {
262 // vector<shape_t*> cryst = load_shapes("/ecl/data/crystal_shape_forward.dat");
263 vector<shape_t*> cryst = load_shapes(m_sap, ECLParts::forward);
264 vector<G4LogicalVolume*> wrapped_crystals;
265 for (auto it = cryst.begin(); it != cryst.end(); it++) {
266 shape_t* s = *it;
267 wrapped_crystals.push_back(wrapped_crystal(s, "forward", 0.20 - 0.02));
268 }
269
270 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
271 const cplacement_t& t = *it;
272 auto s = find_if(cryst.begin(), cryst.end(), [&t](const shape_t* shape) {return shape->nshape == t.nshape;});
273 if (s == cryst.end()) continue;
274
275 G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(t);
276 int indx = it - bp.begin();
277 new G4PVPlacement(twc, wrapped_crystals[s - cryst.begin()], suf("ECLForwardWrappedCrystal_Physical", indx), crystalvolume_logical,
278 false, indx, overlap);
279 }
280
281 for (shape_t* shape : cryst) {
282 delete shape;
283 }
284 }
285
286 if (b_preamplifier) {
287 for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
288 G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(*it);
289 int indx = it - bp.begin();
290 auto pv = new G4PVPlacement(twc * G4TranslateZ3D(300 / 2 + 0.20 + get_pa_box_height() / 2)*G4RotateZ3D(-M_PI / 2), get_preamp(),
291 suf("phys_forward_preamplifier", indx), crystalvolume_logical, false, indx, 0);
292 if (overlap)pv->CheckOverlaps(1000);
293 }
294 }
295
296 if (b_support_leg) {
297 const G4VisAttributes* batt = att("iron");
298
299 G4VSolid* s1 = new G4Box("fwd_leg_p1", 130. / 2, 170. / 2, (40. - milled_thickness) / 2);
300 G4LogicalVolume* l1 = new G4LogicalVolume(s1, Materials::get("SUS304"), "l1", 0, 0, 0);
301 G4Transform3D t1 = G4Translate3D(0, 170. / 2, (40. - milled_thickness) / 2 + milled_thickness);
302 l1->SetVisAttributes(batt);
303
304 G4VSolid* s2 = new G4Box("fwd_leg_p2", 60. / 2, 130. / 2, 137. / 2);
305 G4LogicalVolume* l2 = new G4LogicalVolume(s2, Materials::get("SUS304"), "l2", 0, 0, 0);
306 G4Transform3D t2 = G4Translate3D(0, 130. / 2 + 35, 40. + 137. / 2);
307 l2->SetVisAttributes(batt);
308
309 Point_t v3[] = {{ -75. / 2, -265. / 2}, {75. / 2 - 30, -265. / 2}, {75. / 2, -265. / 2 + 30}, {75. / 2, 265. / 2}, { -75. / 2, 265. / 2}};
310 const int n3 = sizeof(v3) / sizeof(Point_t);
311 G4ThreeVector c3[n3 * 2];
312
313 for (int i = 0; i < n3; i++) c3[i + 0] = G4ThreeVector(v3[i].x, v3[i].y, -140. / 2);
314 for (int i = 0; i < n3; i++) c3[i + n3] = G4ThreeVector(v3[i].x, v3[i].y, 140. / 2);
315
316 G4VSolid* s3 = new BelleCrystal("fwd_leg_p3", n3, c3);
317 G4LogicalVolume* l3 = new G4LogicalVolume(s3, Materials::get("SUS304"), "l3", 0, 0, 0);
318 G4Transform3D t3 = G4Translate3D(0, 265. / 2 + 35, 40. + 137. + 75. / 2) * G4RotateY3D(-M_PI / 2);
319 l3->SetVisAttributes(batt);
320
321 G4VSolid* s4 = new G4Box("fwd_leg_p4", 130. / 2, 5. / 2, (5. + milled_thickness) / 2);
322 G4LogicalVolume* l4 = new G4LogicalVolume(s4, Materials::get("SUS304"), "l4", 0, 0, 0);
323 G4Transform3D t4 = G4Translate3D(0, 170. - 5. / 2, (milled_thickness - 5.) / 2);
324 l4->SetVisAttributes(batt);
325
326 G4VSolid* s5 = new G4Box("fwd_leg_p5", 140. / 2, 130. / 2, 80. / 2);
327 G4LogicalVolume* l5 = new G4LogicalVolume(s5, Materials::get("SUS304"), "l5", 0, 0, 0);
328 G4Transform3D t5 = G4Translate3D(0, 180. + 130. / 2, 97. + 80. / 2);
329 l5->SetVisAttributes(batt);
330
331 G4VSolid* s6 = new G4Box("fwd_leg_p6", 140. / 2, 110. / 2, 160. / 2);
332 G4LogicalVolume* l6 = new G4LogicalVolume(s6, Materials::get("G4_AIR"), "l6", 0, 0, 0);
333 G4Transform3D t6 = G4Translate3D(0, 310. + 110. / 2, 97. + 160. / 2);
334 l6->SetVisAttributes(att("air"));
335
336 G4VSolid* s6a = new G4Box("fwd_leg_p6a", 140. / 2, (110. - 45.) / 2, 160. / 2);
337 G4LogicalVolume* l6a = new G4LogicalVolume(s6a, Materials::get("SUS304"), "l6a", 0, 0, 0);
338 l6a->SetVisAttributes(batt);
339 new G4PVPlacement(G4TranslateY3D(-45. / 2), l6a, "l6a_physical", l6, false, 0, overlap);
340
341 G4VSolid* s6b = new G4Box("fwd_leg_p6b", 60. / 2, 45. / 2, 160. / 2);
342 G4LogicalVolume* l6b = new G4LogicalVolume(s6b, Materials::get("SUS304"), "l6b", 0, 0, 0);
343 l6b->SetVisAttributes(batt);
344 double dy = 110. / 2 - 45 + 45. / 2;
345 new G4PVPlacement(G4TranslateY3D(dy), l6b, "l6b_physical", l6, false, 0, overlap);
346
347 G4VSolid* s6c = new G4Box("fwd_leg_p6c", 40. / 2, 45. / 2, 22.5 / 2);
348 G4LogicalVolume* l6c = new G4LogicalVolume(s6c, Materials::get("SUS304"), "l6c", 0, 0, 0);
349 l6c->SetVisAttributes(batt);
350 new G4PVPlacement(G4Translate3D(30 + 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 0, overlap);
351 new G4PVPlacement(G4Translate3D(30 + 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 1, overlap);
352 new G4PVPlacement(G4Translate3D(-30 - 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 2, overlap);
353 new G4PVPlacement(G4Translate3D(-30 - 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 3, overlap);
354
355 G4AssemblyVolume* support_leg = new G4AssemblyVolume();
356
357 support_leg->AddPlacedVolume(l1, t1);
358 support_leg->AddPlacedVolume(l2, t2);
359 support_leg->AddPlacedVolume(l3, t3);
360 support_leg->AddPlacedVolume(l4, t4);
361 support_leg->AddPlacedVolume(l5, t5);
362 support_leg->AddPlacedVolume(l6, t6);
363
364 for (int i = 0; i < 8; i++) {
365 G4Transform3D tp = gTrans * G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4) *
366 G4Translate3D(0, 1415 - 165 + 420. / 2, ZT + (97. + 160.) / 2) * G4Translate3D(0, -420. / 2, -(97. + 160.) / 2);
367 support_leg->MakeImprint(top, tp, 0, overlap);
368 }
369
370 // G4VSolid* s_all = new G4Box("leg_all", 140. / 2, 420. / 2, (97. + 160) / 2);
371 // G4LogicalVolume* l_all = new G4LogicalVolume(s_all, Materials::get("G4_AIR"), "l_all", 0, 0, 0);
372 // l_all->SetVisAttributes(att("silv"));
373 // G4Transform3D tp = G4Translate3D(0, -420. / 2, -(97. + 160.) / 2);
374 // support_leg->MakeImprint(l_all, tp, 0, overlap);
375
376 // for (int i = 0; i < 8; i++)
377 // new G4PVPlacement(G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0, 1415 - 165 + 420. / 2,
378 // 1960 + ZT + (97. + 160.) / 2), l_all, suf("support_leg_physical", i), top, false, i, overlap);
379
380 // for (int i = 0; i < 8; i++)
381 // new G4PVPlacement(G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0, 1415 - 165 + 420. / 2,
382 // 1960 + ZT + (97. + 160.) / 2)*tp * t4, l4, suf("support_leg_p4_physical", i), top, false, i, overlap);
383
384 delete support_leg;
385 }
386
387
388 if (b_support_structure_13) { // numbering scheme as in ECL-004K102.pdf page 13
389
390 // Define one layer as one assembly volume
391 G4AssemblyVolume* acs = new G4AssemblyVolume();
392
393 double Z0 = 434, Ro = 1415 - 20;
394
395 // double R2_0 = 667, R2_1 = 1245/cos(M_PI/16), dR2 = R2_1 - R2_0, dz = dR2*cos(M_PI/16);
396 // G4VSolid* sv_crystal_support1 = new G4Trd("sv_crystal_support1", 30/2, 30/2, R2_0*sin(M_PI/16)-40*cos(M_PI/16), R2_1*sin(M_PI/16)-40*cos(M_PI/16), dz/2);
397 // G4LogicalVolume* lv_crystal_support1 = new G4LogicalVolume(sv_crystal_support1, world_mat, "lv_crystal_support1", 0, 0, 0);
398 // // lv_crystal_support->SetVisAttributes(airvol);
399 // new G4PVPlacement(G4Translate3D(R2_0*cos(M_PI/16)+dz/2+3,0,Z0-95)*G4RotateY3D(M_PI/2), lv_crystal_support1, "phys_crystal_support1", crystalSectorLogical, false, 0, overlaps);
400
401 // double dx=425+33, dy1 = dx*tan(M_PI/16), dz = 630+20-67, dy2 = dy1+dz*tan(M_PI/16);
402 // G4VSolid* sv_crystal_support1 = new G4Trd("sv_crystal_support1", 30/2, 30/2, dy1, dy2, dz/2);
403 // G4LogicalVolume* lv_crystal_support1 = new G4LogicalVolume(sv_crystal_support1, world_mat, "lv_crystal_support1", 0, 0, 0);
404 // // lv_crystal_support->SetVisAttributes(airvol);
405 // new G4PVPlacement(G4Translate3D(40/sin(M_PI/16)+dx+dz/2,0,Z0-95)*G4RotateY3D(M_PI/2), lv_crystal_support1, "phys_crystal_support1", crystalSectorLogical, false, 0, overlaps);
406
407 G4VSolid* solid10_p1 = new G4Box("fwd_solid10_p1", 558. / 2 + 9.5, 20. / 2, 105. / 2);
408 G4VSolid* solid10_p2 = new G4Box("fwd_solid10_p2", 559. / 2 + 9.5, 11. / 2, 71. / 2);
409 G4VSolid* solid10_p3 = new G4SubtractionSolid("fwd_solid10_p3", solid10_p1, solid10_p2, G4Translate3D(0, -11. / 2, 71. / 2 - 17.5));
410
411 G4LogicalVolume* lsolid10 = new G4LogicalVolume(solid10_p3, Materials::get("A6063"), "lsolid10", 0, 0, 0);
412 lsolid10->SetVisAttributes(att("alum"));
413 G4Transform3D tsolid10_p1(G4RotateZ3D(M_PI / 16)*G4Translate3D(954.5 + 2.55 - 1, -30, Z0 - 105. / 2 - 5));
414 acs->AddPlacedVolume(lsolid10, tsolid10_p1);
415 // new G4PVPlacement(tsolid10_p1, lsolid10, "psolid10_p1", crystalSectorLogical, false, 0, overlaps);
416 G4Transform3D tsolid10_p2(G4RotateZ3D(-M_PI / 16)*G4Translate3D(954.5 + 2.55 - 1, 30, Z0 - 105. / 2 - 5)*G4RotateZ3D(M_PI));
417 acs->AddPlacedVolume(lsolid10, tsolid10_p2);
418 // new G4PVPlacement(tsolid10_p2, lsolid10, "psolid10_p2", crystalSectorLogical, false, 0, overlaps);
419
420 G4VSolid* solid1_p1 = new G4Box("fwd_solid1_p1", 100. / 2, 30. / 2, 30. / 2);
421 G4VSolid* solid1_p2 = new G4Box("fwd_solid1_p2", 80. / 2, 10. / 2, 31. / 2);
422 G4VSolid* solid1_p3 = new G4SubtractionSolid("fwd_solid1_p3", solid1_p1, solid1_p2, G4Transform3D::Identity);
423
424 G4LogicalVolume* lsolid1 = new G4LogicalVolume(solid1_p3, Materials::get("A5052"), "lsolid1", 0, 0, 0);
425 lsolid1->SetVisAttributes(att("alum"));
426
427 G4Transform3D tsolid1_p1(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 3 * 140, 0, Z0 - 95));
428 acs->AddPlacedVolume(lsolid1, tsolid1_p1);
429 // new G4PVPlacement(tsolid1_p1, lsolid1, "psolid1_p1", crystalSectorLogical, false, 0, overlaps);
430 G4Transform3D tsolid1_p2(G4RotateZ3D(M_PI / 16 * (1 - 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 3 * 140, 0, Z0 - 95));
431 acs->AddPlacedVolume(lsolid1, tsolid1_p2);
432 // new G4PVPlacement(tsolid1_p2, lsolid1, "psolid1_p2", crystalSectorLogical, false, 0, overlaps);
433
434 G4Transform3D tsolid1_p3(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 2 * 140, 0, Z0 - 95));
435 acs->AddPlacedVolume(lsolid1, tsolid1_p3);
436 // new G4PVPlacement(tsolid1_p3, lsolid1, "psolid1_p3", crystalSectorLogical, false, 0, overlaps);
437 G4Transform3D tsolid1_p4(G4RotateZ3D(M_PI / 16 * (1 - 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 2 * 140, 0, Z0 - 95));
438 acs->AddPlacedVolume(lsolid1, tsolid1_p4);
439 // new G4PVPlacement(tsolid1_p4, lsolid1, "psolid1_p4", crystalSectorLogical, false, 0, overlaps);
440
441 G4Transform3D tsolid1_p5(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 1 * 140, 0, Z0 - 95));
442 acs->AddPlacedVolume(lsolid1, tsolid1_p5);
443 // new G4PVPlacement(tsolid1_p5, lsolid1, "psolid1_p5", crystalSectorLogical, false, 0, overlaps);
444 G4Transform3D tsolid1_p6(G4RotateZ3D(M_PI / 16 * (1 - 2 / 3.))*G4Translate3D(Ro - 8 - 50 - 1 * 140, 0, Z0 - 95));
445 acs->AddPlacedVolume(lsolid1, tsolid1_p6);
446 // new G4PVPlacement(tsolid1_p6, lsolid1, "psolid1_p6", crystalSectorLogical, false, 0, overlaps);
447
448 G4Transform3D tsolid1_p8(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3.))*G4Translate3D(Ro - 8 - 50, 0, Z0 - 100));
449 acs->AddPlacedVolume(lsolid1, tsolid1_p8);
450 // new G4PVPlacement(tsolid1_p8, lsolid1, "psolid1_p8", crystalSectorLogical, false, 0, overlaps);
451 G4Transform3D tsolid1_p9(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3. + (1 + 1. / 3) / 3))*G4Translate3D(Ro - 8 - 50, 0, Z0 - 100));
452 acs->AddPlacedVolume(lsolid1, tsolid1_p9);
453 // new G4PVPlacement(tsolid1_p9, lsolid1, "psolid1_p9", crystalSectorLogical, false, 0, overlaps);
454 G4Transform3D tsolid1_p10(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3. + 2 * (1 + 1. / 3) / 3))*G4Translate3D(Ro - 8 - 50, 0, Z0 - 100));
455 acs->AddPlacedVolume(lsolid1, tsolid1_p10);
456 // new G4PVPlacement(tsolid1_p10, lsolid1, "psolid1_p10", crystalSectorLogical, false, 0, overlaps);
457
458 G4VSolid* solid1_p11 = new G4Box("fwd_solid1_p1", 100. / 2, 10. / 2, 30. / 2);
459 G4LogicalVolume* lsolid1_p2 = new G4LogicalVolume(solid1_p11, Materials::get("A5052"), "lsolid1_p2", 0, 0, 0);
460 lsolid1_p2->SetVisAttributes(att("alum"));
461 G4Transform3D tsolid1_p11(G4RotateZ3D(M_PI / 16 * (-1 + 2 / 3. - 1. / 3))*G4Translate3D(Ro - 8 - 50, 0, Z0 - 100));
462 acs->AddPlacedVolume(lsolid1_p2, tsolid1_p11);
463 // new G4PVPlacement(tsolid1_p11, lsolid1_p2, "psolid1_p11", crystalSectorLogical, false, 0, overlaps);
464
465 G4VSolid* solid1_p12 = new G4Box("fwd_solid1_p1", 86. / 2, 10. / 2, 30. / 2);
466 G4LogicalVolume* lsolid1_p3 = new G4LogicalVolume(solid1_p12, Materials::get("A5052"), "lsolid1_p3", 0, 0, 0);
467 lsolid1_p3->SetVisAttributes(att("alum"));
468 double alpha_p12 = M_PI / 16 * (-1 + 1. / 3);
469 G4Transform3D tsolid1_p12(G4Translate3D(532.2 + 43, 0, Z0 - 75));
470 acs->AddPlacedVolume(lsolid1_p3, tsolid1_p12);
471 // new G4PVPlacement(tsolid1_p12, lsolid1_p3, "psolid1_p12", crystalSectorLogical, false, 0, overlaps);
472 G4Transform3D tsolid1_p13(G4RotateZ3D(alpha_p12)*G4Translate3D(532.2 + 43, 0, Z0 - 75));
473 acs->AddPlacedVolume(lsolid1_p3, tsolid1_p13);
474 // new G4PVPlacement(tsolid1_p13, lsolid1_p3, "psolid1_p13", crystalSectorLogical, false, 0, overlaps);
475 G4Transform3D tsolid1_p14(G4RotateZ3D(-alpha_p12)*G4Translate3D(532.2 + 43, 0, Z0 - 75));
476 acs->AddPlacedVolume(lsolid1_p3, tsolid1_p14);
477 // new G4PVPlacement(tsolid1_p14, lsolid1_p3, "psolid1_p14", crystalSectorLogical, false, 0, overlaps);
478
479 G4VSolid* solid1_p4 = new G4Box("fwd_solid1_p4", 160. / 2, 30. / 2, 30. / 2);
480 G4VSolid* solid1_p5 = new G4Box("fwd_solid1_p5", 140. / 2, 10. / 2, 31. / 2);
481 G4VSolid* solid1_p7 = new G4SubtractionSolid("fwd_solid1_p7", solid1_p4, solid1_p5, G4Transform3D::Identity);
482 G4LogicalVolume* lsolid1_p7 = new G4LogicalVolume(solid1_p7, Materials::get("A5052"), "lsolid1_p7", 0, 0, 0);
483 lsolid1_p7->SetVisAttributes(att("alum"));
484 G4Transform3D tsolid1_p7(G4Translate3D(Ro - 8 - 80 - 4 * 140 + 4, 0, Z0 - 95));
485 acs->AddPlacedVolume(lsolid1_p7, tsolid1_p7);
486 // new G4PVPlacement(tsolid1_p7, lsolid1_p7, "psolid1_p7", crystalSectorLogical, false, 0, overlaps);
487
488 // int sol2count = 0;
489 auto get_bracket = [&](double L, double ang, G4Transform3D & lt) {
490 double thick = 3;
491 double dL = (ang > 0) ? 0 : thick * abs(tan(ang));
492
493 G4VSolid* solid2_p1 = new G4Box("fwd_solid2_p1", (L - 2 * dL) / 2, thick / 2, 30. / 2);
494 G4VSolid* solid2_p2 = new G4Box("fwd_solid2_p2", thick / 2, (15. - dL) / 2, 30. / 2);
495 double dx = thick / 2, y0 = (15. + dL) / 2;
496 G4Transform3D t1(G4Translate3D(L / 2, -dx, 0.)*G4RotateZ3D(-ang)*G4Translate3D(-dx, y0, 0));
497 G4Transform3D t2(G4Translate3D(-L / 2, -dx, 0.)*G4RotateZ3D(ang)*G4Translate3D(dx, y0, 0));
498 G4VSolid* solid2_p3 = new G4UnionSolid("fwd_solid2_p3", solid2_p1, solid2_p2, t1);
499 G4VSolid* solid2_p4 = new G4UnionSolid("fwd_solid2_p4", solid2_p3, solid2_p2, t2);
500
501 G4Transform3D u((ang > 0) ? G4Transform3D::Identity : G4RotateZ3D(M_PI));
502 lt = u * G4Translate3D(dx, 0, 0) * G4RotateZ3D(-M_PI / 2);
503
504 return solid2_p4;
505 };
506 double obj2_dz = Z0 - 95;
507 auto place_solid2 = [&](double dz, double L, double ang, double phi, double mx, double dy) {
508 G4Transform3D lt;
509 G4LogicalVolume* lsolid2 = new G4LogicalVolume(get_bracket(L, ang, lt), Materials::get("A5052"), "lsolid2", 0, 0, 0);
510 lsolid2->SetVisAttributes(att("alum"));
511
512 G4Transform3D tsolid2_p1(G4RotateZ3D(phi)*G4Translate3D(mx, dy, dz)*lt);
513 // string pname("psolid2_p"); pname += to_string(++sol2count);
514 acs->AddPlacedVolume(lsolid2, tsolid2_p1);
515 // new G4PVPlacement(tsolid2_p1, lsolid2, pname.c_str(), crystalSectorLogical, false, 0, overlaps);
516 };
517 auto place_solid3 = [&](double L, double ang, const G4Transform3D & t) {
518 G4Transform3D lt;
519 G4LogicalVolume* lsolid2 = new G4LogicalVolume(get_bracket(L, ang, lt), Materials::get("A5052"), "lsolid2", 0, 0, 0);
520 lsolid2->SetVisAttributes(att("alum"));
521
522 G4Transform3D tsolid2_p1(t * lt);
523 // string pname("psolid2_p"); pname += to_string(++sol2count);
524 acs->AddPlacedVolume(lsolid2, tsolid2_p1);
525 // new G4PVPlacement(tsolid2_p1, lsolid2, pname.c_str(), crystalSectorLogical, false, 0, overlaps);
526 };
527
528 G4Point3D aa(-50 + 15, 15, 0), bb(-50 + 15, -15, 0);
529
530 aa = tsolid1_p1 * G4Point3D(-50 + 15, 15, 0);
531 bb = tsolid1_p2 * G4Point3D(-50 + 15, -15, 0);
532 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
533
534 aa = tsolid1_p1 * G4Point3D(50 - 15, 15, 0);
535 bb = tsolid1_p2 * G4Point3D(50 - 15, -15, 0);
536 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
537
538 G4Point3D r0(0, 40 / cos(M_PI / 16), 0);
539 G4Vector3D np(sin(M_PI / 16), cos(M_PI / 16), 0), mid;
540 G4RotateZ3D rb(M_PI / 24);
541 double L, phi_uu;
542
543
544 phi_uu = (tsolid1_p1 * G4Vector3D(1, 0, 0)).angle(G4Vector3D(cos(M_PI / 16), -sin(M_PI / 16), 0));
545 G4Vector3D n = tsolid1_p1 * G4Vector3D(-sin(phi_uu / 2), -cos(phi_uu / 2), 0);
546 double xj = 50 - 15;
547 aa = tsolid1_p1 * G4Point3D(-xj, -15, 0);
548 L = -((aa - r0) * np) / (n * np);
549 place_solid3(L, -phi_uu / 2, tsolid1_p1 * G4Translate3D(-xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
550 place_solid3(L, -phi_uu / 2, tsolid1_p2 * G4Translate3D(-xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
551 aa = tsolid1_p1 * G4Point3D(xj, -15, 0);
552 L = -((aa - r0) * np) / (n * np);
553 place_solid3(L, phi_uu / 2, tsolid1_p1 * G4Translate3D(xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
554 place_solid3(L, phi_uu / 2, tsolid1_p2 * G4Translate3D(xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
555
556 aa = tsolid1_p3 * G4Point3D(-xj, -15, 0);
557 L = -((aa - r0) * np) / (n * np);
558 place_solid3(L, -phi_uu / 2, tsolid1_p3 * G4Translate3D(-xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
559 place_solid3(L, -phi_uu / 2, tsolid1_p4 * G4Translate3D(-xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
560 aa = tsolid1_p3 * G4Point3D(xj, -15, 0);
561 L = -((aa - r0) * np) / (n * np);
562 place_solid3(L, phi_uu / 2, tsolid1_p3 * G4Translate3D(xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
563 place_solid3(L, phi_uu / 2, tsolid1_p4 * G4Translate3D(xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
564
565 aa = tsolid1_p5 * G4Point3D(-xj, -15, 0);
566 L = -((aa - r0) * np) / (n * np);
567 place_solid3(L, -phi_uu / 2, tsolid1_p5 * G4Translate3D(-xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
568 place_solid3(L, -phi_uu / 2, tsolid1_p6 * G4Translate3D(-xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
569 aa = tsolid1_p5 * G4Point3D(xj, -15, 0);
570 L = -((aa - r0) * np) / (n * np);
571 place_solid3(L, phi_uu / 2, tsolid1_p5 * G4Translate3D(xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
572 place_solid3(L, phi_uu / 2, tsolid1_p6 * G4Translate3D(xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
573
574 phi_uu = M_PI / 16;
575 n = tsolid1_p7 * G4Vector3D(-sin(phi_uu / 2), -cos(phi_uu / 2), 0);
576 xj = 80 - 15;
577 aa = tsolid1_p7 * G4Point3D(-xj, -15, 0);
578 L = -((aa - r0) * np) / (n * np);
579 place_solid3(L, -phi_uu / 2, tsolid1_p7 * G4Translate3D(-xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
580 place_solid3(L, -phi_uu / 2, tsolid1_p7 * G4Translate3D(-xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
581 aa = tsolid1_p7 * G4Point3D(xj, -15, 0);
582 L = -((aa - r0) * np) / (n * np);
583 place_solid3(L, phi_uu / 2, tsolid1_p7 * G4Translate3D(xj, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
584 place_solid3(L, phi_uu / 2, tsolid1_p7 * G4Translate3D(xj, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
585
586
587 aa = tsolid1_p3 * G4Point3D(-50 + 15, 15, 0);
588 bb = tsolid1_p4 * G4Point3D(-50 + 15, -15, 0);
589 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
590
591 aa = tsolid1_p3 * G4Point3D(50 - 15, 15, 0);
592 bb = tsolid1_p4 * G4Point3D(50 - 15, -15, 0);
593 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
594
595
596 aa = tsolid1_p5 * G4Point3D(-50 + 15, 15, 0);
597 bb = tsolid1_p6 * G4Point3D(-50 + 15, -15, 0);
598 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
599
600 aa = tsolid1_p5 * G4Point3D(50 - 15, 15, 0);
601 bb = tsolid1_p6 * G4Point3D(50 - 15, -15, 0);
602 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
603
604
605 obj2_dz = Z0 - 100;
606 aa = tsolid1_p8 * G4Point3D(-50 + 15, 15, 0);
607 bb = tsolid1_p9 * G4Point3D(-50 + 15, -15, 0);
608 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 72, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
609
610 aa = tsolid1_p8 * G4Point3D(50 - 15, 15, 0);
611 bb = tsolid1_p9 * G4Point3D(50 - 15, -15, 0);
612 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 72, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
613
614 aa = tsolid1_p9 * G4Point3D(-50 + 15, 15, 0);
615 bb = tsolid1_p10 * G4Point3D(-50 + 15, -15, 0);
616 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 72, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
617
618 aa = tsolid1_p9 * G4Point3D(50 - 15, 15, 0);
619 bb = tsolid1_p10 * G4Point3D(50 - 15, -15, 0);
620 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 72, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
621
622 phi_uu = (tsolid1_p8 * G4Vector3D(1, 0, 0)).angle(tsolid1_p11 * G4Vector3D(1, 0, 0));
623 r0 = tsolid1_p11 * G4Point3D(-50 + 15, 5, 0);
624 np = tsolid1_p11 * G4Vector3D(0, 1, 0);
625
626 n = tsolid1_p8 * G4Vector3D(-sin(phi_uu / 2), -cos(phi_uu / 2), 0);
627 aa = tsolid1_p8 * G4Point3D(-50 + 15, -15, 0);
628 L = -((aa - r0) * np) / (n * np);
629 place_solid3(L, -phi_uu / 2, tsolid1_p8 * G4Translate3D(-50 + 15, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
630
631 aa = tsolid1_p8 * G4Point3D(50 - 15, -15, 0);
632 L = -((aa - r0) * np) / (n * np);
633 place_solid3(L, phi_uu / 2, tsolid1_p8 * G4Translate3D(50 - 15, -15, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
634
635 phi_uu = (tsolid1_p11 * G4Vector3D(1, 0, 0)).angle(G4RotateZ3D(-M_PI / 16) * G4Vector3D(1, 0, 0));
636 r0 = G4RotateZ3D(-M_PI / 16) * G4Point3D(0, 40, 0);
637 np = G4RotateZ3D(-M_PI / 16) * G4Vector3D(0, 1, 0);
638
639 n = tsolid1_p11 * G4Vector3D(-sin(phi_uu / 2), -cos(phi_uu / 2), 0);
640 aa = tsolid1_p11 * G4Point3D(-50 + 15, -5, 0);
641 L = -((aa - r0) * np) / (n * np);
642 place_solid3(L, -phi_uu / 2, tsolid1_p11 * G4Translate3D(-50 + 15, -5, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
643
644 aa = tsolid1_p11 * G4Point3D(50 - 15, -5, 0);
645 L = -((aa - r0) * np) / (n * np);
646 place_solid3(L, phi_uu / 2, tsolid1_p11 * G4Translate3D(50 - 15, -5, 0)*G4RotateZ3D(-phi_uu / 2)*G4Translate3D(0, -L / 2, 0));
647
648 G4Transform3D ttt(G4RotateZ3D(M_PI / 16)*G4Translate3D(Ro - 8 - 50, 0, Z0 - 100));
649 phi_uu = (tsolid1_p10 * G4Vector3D(1, 0, 0)).angle(ttt * G4Vector3D(1, 0, 0));
650 r0 = ttt * G4Point3D(0, -40, 0);
651 np = ttt * G4Vector3D(0, 1, 0);
652 n = tsolid1_p10 * G4Vector3D(-sin(phi_uu / 2), cos(phi_uu / 2), 0);
653 aa = tsolid1_p10 * G4Point3D(-50 + 15, 15, 0);
654 L = -((aa - r0) * np) / (n * np);
655 place_solid3(L, -phi_uu / 2, tsolid1_p10 * G4Translate3D(-50 + 15, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
656
657 aa = tsolid1_p10 * G4Point3D(50 - 15, 15, 0);
658 L = -((aa - r0) * np) / (n * np);
659 place_solid3(L, phi_uu / 2, tsolid1_p10 * G4Translate3D(50 - 15, 15, 0)*G4RotateZ3D(phi_uu / 2)*G4Translate3D(0, L / 2, 0));
660
661 obj2_dz = Z0 - 75;
662 aa = tsolid1_p12 * G4Point3D(-43 + 15, -5, 0);
663 bb = tsolid1_p13 * G4Point3D(-43 + 15, 5, 0);
664 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
665
666 aa = tsolid1_p12 * G4Point3D(43 - 15, -5, 0);
667 bb = tsolid1_p13 * G4Point3D(43 - 15, 5, 0);
668 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
669
670 aa = tsolid1_p12 * G4Point3D(-43 + 15, 5, 0);
671 bb = tsolid1_p14 * G4Point3D(-43 + 15, -5, 0);
672 place_solid2(obj2_dz, (aa - bb).mag(), -M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
673
674 aa = tsolid1_p12 * G4Point3D(43 - 15, 5, 0);
675 bb = tsolid1_p14 * G4Point3D(43 - 15, -5, 0);
676 place_solid2(obj2_dz, (aa - bb).mag(), M_PI / 48, (aa + bb).phi(), (aa + bb).rho() / 2, 0);
677
678 G4VSolid* solid11_p1 = new G4Box("fwd_solid11_p1", 80. / 2, 30. / 2, 40. / 2);
679 G4VSolid* solid11_p2 = new G4Box("fwd_solid11_p2", 81. / 2, 30. / 2, 40. / 2);
680 G4VSolid* solid11_p3 = new G4SubtractionSolid("fwd_solid11_p3", solid11_p1, solid11_p2, G4Translate3D(0, -3, 3));
681
682 G4LogicalVolume* lsolid11 = new G4LogicalVolume(solid11_p3, Materials::get("SUS304"), "lsolid11", 0, 0, 0);
683 lsolid11->SetVisAttributes(att("iron"));
684 G4Transform3D tsolid11_p1(G4RotateZ3D(M_PI / 16)*G4Translate3D(580, -35, Z0 - 40));
685 acs->AddPlacedVolume(lsolid11, tsolid11_p1);
686
687 G4Transform3D tsolid11_p2(G4RotateZ3D(-M_PI / 16)*G4Translate3D(580, 35, Z0 - 40)*G4RotateZ3D(M_PI));
688 acs->AddPlacedVolume(lsolid11, tsolid11_p2);
689
690 G4VSolid* solid12_p1 = new G4Box("fwd_solid12_p1", 120. / 2, 20. / 2, 115. / 2);
691 G4VSolid* solid12_p2 = new G4Box("fwd_solid12_p2", 121. / 2, 12. / 2, 86. / 2);
692 G4VSolid* solid12_p3 = new G4SubtractionSolid("fwd_solid12_p3", solid12_p1, solid12_p2, G4Translate3D(0, -11. / 2 + 1,
693 86. / 2 - 17.5 - 5));
694 G4LogicalVolume* lsolid12 = new G4LogicalVolume(solid12_p3, Materials::get("A6063"), "lsolid12", 0, 0, 0);
695 lsolid12->SetVisAttributes(att("alum"));
696
697 G4Transform3D tsolid12_p2(G4RotateZ3D(-M_PI / 16)*G4Translate3D(Ro - 8 - 60, 30, Z0 - 115. / 2)*G4RotateZ3D(M_PI));
698 acs->AddPlacedVolume(lsolid12, tsolid12_p2);
699
700 G4VSolid* solid12r_p1 = new G4Box("fwd_solid12r_p1", 100. / 2, 30. / 2, 75. / 2);
701 G4VSolid* solid12r_p2 = new G4Box("fwd_solid12r_p2", 101. / 2, 21. / 2, 41. / 2);
702 G4VSolid* solid12r_p3 = new G4SubtractionSolid("fwd_solid12r_p3", solid12r_p1, solid12r_p2, G4Translate3D(0, -21. / 2 + 5,
703 -41. / 2 - 75. / 2 + 40));
704 G4LogicalVolume* lsolid12r = new G4LogicalVolume(solid12r_p3, Materials::get("A6063"), "lsolid12r", 0, 0, 0);
705 lsolid12r->SetVisAttributes(att("alum"));
706 G4Transform3D tsolid12r_p1(G4RotateZ3D(M_PI / 16)*G4Translate3D(Ro - 8 - 50, -30 - 15, Z0 - 40 - 75. / 2));
707 acs->AddPlacedVolume(lsolid12r, tsolid12r_p1);
708
709 G4Transform3D tr = G4RotateZ3D(M_PI / 16) * G4ReflectY3D();
710 acs->MakeImprint(innervolumesector_logical, tr, 0, overlap);
711 tr = G4RotateZ3D(-M_PI / 16);
712 acs->MakeImprint(innervolumesector_logical, tr, 1, overlap);
713
714 delete acs;
715 }
716
717 if (b_support_structure_15) { // numbering scheme as in ECL-004K102.pdf page 15
718 G4AssemblyVolume* acs = new G4AssemblyVolume();
719
720 double Z0 = 434, Ro = 1415 - 20;
721
722 G4VSolid* solid1_p1 = new G4Box("fwd_solid1_p1", 10. / 2, 398. / 2, 5. / 2);
723 G4VSolid* solid1_p2 = new G4Box("fwd_solid1_p2", 4. / 2, 398. / 2 - 32, 6. / 2);
724 G4VSolid* solid1_p3 = new G4SubtractionSolid("fwd_solid1_p3", solid1_p1, solid1_p2, G4Transform3D(G4RotationMatrix(),
725 G4ThreeVector(-4,
726 0, 0)));
727
728 G4LogicalVolume* lsolid1 = new G4LogicalVolume(solid1_p3, Materials::get("SUS304"), "lsolid1", 0, 0, 0);
729 lsolid1->SetVisAttributes(att("iron"));
730 G4Transform3D tsolid1_p1(G4Translate3D(1350, -16, Z0 - 40 + 3 + 2.5));
731 acs->AddPlacedVolume(lsolid1, tsolid1_p1);
732
733 G4VSolid* solid1a_p1 = new G4Box("fwd_solid1a_p1", 10. / 2, 348.5 / 2, 5. / 2);
734 G4VSolid* solid1a_p2 = new G4Box("fwd_solid1a_p2", 4. / 2, 348.5 / 2 - 32, 6. / 2);
735 G4VSolid* solid1a_p3 = new G4SubtractionSolid("fwd_solid1a_p3", solid1a_p1, solid1a_p2, G4Transform3D(G4RotationMatrix(),
736 G4ThreeVector(4, 0, 0)));
737
738 G4LogicalVolume* lsolid1a = new G4LogicalVolume(solid1a_p3, Materials::get("SUS304"), "lsolid1a", 0, 0, 0);
739 lsolid1a->SetVisAttributes(att("iron"));
740 G4Transform3D tsolid1a_p1(G4Translate3D(1250, -16, Z0 - 40 + 3 + 2.5));
741 acs->AddPlacedVolume(lsolid1a, tsolid1a_p1);
742
743 G4VSolid* solid1b_p1 = new G4Box("fwd_solid1b_p1", 10. / 2, (210 + 210 + 16) / 2, 5. / 2);
744 G4VSolid* solid1b_p2 = new G4Box("fwd_solid1b_p2", 4. / 2, (210 + 210 + 16) / 2 - 16, 6. / 2);
745 G4VSolid* solid1b_p3 = new G4SubtractionSolid("fwd_solid1b_p3", solid1b_p1, solid1b_p2, G4Transform3D(G4RotationMatrix(),
746 G4ThreeVector(4, 0, 0)));
747
748 G4LogicalVolume* lsolid1b = new G4LogicalVolume(solid1b_p3, Materials::get("SUS304"), "lsolid1b", 0, 0, 0);
749 lsolid1b->SetVisAttributes(att("iron"));
750 G4Transform3D t1b(G4Translate3D(1204 - (210 + 210 + 16) / 2 + 16, 0, Z0 - 52 + 8 + 5. / 2));
751 G4Transform3D tsolid1b_p1(t1b * G4TranslateY3D(-10)*G4RotateZ3D(-M_PI / 2));
752 acs->AddPlacedVolume(lsolid1b, tsolid1b_p1);
753 G4Transform3D tsolid1b_p2(t1b * G4TranslateY3D(10)*G4RotateZ3D(M_PI / 2));
754 acs->AddPlacedVolume(lsolid1b, tsolid1b_p2);
755 G4Transform3D tsolid1b_p3(t1b * G4TranslateY3D(-110)*G4RotateZ3D(M_PI / 2));
756 acs->AddPlacedVolume(lsolid1b, tsolid1b_p3);
757 G4Transform3D tsolid1b_p4(t1b * G4TranslateY3D(110)*G4RotateZ3D(-M_PI / 2));
758 acs->AddPlacedVolume(lsolid1b, tsolid1b_p4);
759
760 G4VSolid* solid2_p1 = new G4Box("fwd_solid2_p1", 20. / 2, 210. / 2, 5. / 2);
761 G4VSolid* solid2_p2 = new G4Box("fwd_solid2_p2", 14. / 2, 210. / 2 - 16, 6. / 2);
762 G4VSolid* solid2_p3 = new G4SubtractionSolid("fwd_solid2_p3", solid2_p1, solid2_p2, G4Transform3D(G4RotationMatrix(),
763 G4ThreeVector(-4,
764 0, 0)));
765
766 G4LogicalVolume* lsolid2 = new G4LogicalVolume(solid2_p3, Materials::get("SUS304"), "lsolid2", 0, 0, 0);
767 lsolid2->SetVisAttributes(att("iron"));
768 G4Transform3D tsolid2_p1(G4Translate3D(1204 - (210 + 210) - 210 / 2 + 16, 58, Z0 - 52 + 8 + 5. / 2)*G4RotateZ3D(-M_PI / 2));
769 acs->AddPlacedVolume(lsolid2, tsolid2_p1);
770 // new G4PVPlacement(tsolid2_p1, lsolid2, "psolid2_p1", crystalSectorLogical, false, 0, overlap);
771 G4Transform3D tsolid2_p2(G4Translate3D(1204 - (210 + 210) - 210 / 2 + 16, -58, Z0 - 52 + 8 + 5. / 2)*G4RotateZ3D(M_PI / 2));
772 acs->AddPlacedVolume(lsolid2, tsolid2_p2);
773 // new G4PVPlacement(tsolid2_p2, lsolid2, "psolid2_p2", crystalSectorLogical, false, 0, overlap);
774
775 G4VSolid* solid3_p1 = new G4Box("fwd_solid3_p1", 32. / 2, 396. / 2, 8. / 2);
776 G4LogicalVolume* lsolid3 = new G4LogicalVolume(solid3_p1, Materials::get("SUS304"), "lsolid3", 0, 0, 0);
777 lsolid3->SetVisAttributes(att("iron"));
778 G4Transform3D tsolid3_p1(G4Translate3D(1204, 0, Z0 - 52 + 8. / 2));
779 acs->AddPlacedVolume(lsolid3, tsolid3_p1);
780 // new G4PVPlacement(tsolid3_p1, lsolid3, "psolid3_p1", crystalSectorLogical, false, 0, overlap);
781
782 G4VSolid* solid3n_p1 = new G4Box("fwd_solid3n_p1", 32. / 2, 230. / 2, 8. / 2);
783 G4LogicalVolume* lsolid3n = new G4LogicalVolume(solid3n_p1, Materials::get("SUS304"), "lsolid3n", 0, 0, 0);
784 lsolid3n->SetVisAttributes(att("iron"));
785 G4Transform3D tsolid3n_p1(G4Translate3D(1204 - 420, 0, Z0 - 52 + 8. / 2));
786 acs->AddPlacedVolume(lsolid3n, tsolid3n_p1);
787 // new G4PVPlacement(tsolid3n_p1, lsolid3n, "psolid3n_p1", crystalSectorLogical, false, 0, overlap);
788
789 G4VSolid* solid4_p1 = new G4Box("fwd_solid4_p1", 16. / 2, 160. / 2, 8. / 2);
790 G4LogicalVolume* lsolid4 = new G4LogicalVolume(solid4_p1, Materials::get("SUS304"), "lsolid4", 0, 0, 0);
791 lsolid4->SetVisAttributes(att("iron"));
792 G4Transform3D tsolid4_p1(G4Translate3D(598, 0, Z0 - 52 + 8. / 2));
793 acs->AddPlacedVolume(lsolid4, tsolid4_p1);
794 // new G4PVPlacement(tsolid4_p1, lsolid4, "psolid4_p1", crystalSectorLogical, false, 0, overlap);
795
796 G4VSolid* solid5_p1 = new G4Box("fwd_solid5_p1", 650. / 2, 30. / 2, 40. / 2);
797 G4VSolid* solid5_p2 = new G4Box("fwd_solid5_p2", 651. / 2, 30. / 2, 40. / 2);
798 G4VSolid* solid5_p3 = new G4SubtractionSolid("solid5_p3", solid5_p1, solid5_p2, G4Translate3D(0, -3, 3));
799
800 G4LogicalVolume* lsolid5 = new G4LogicalVolume(solid5_p3, Materials::get("SUS304"), "lsolid5", 0, 0, 0);
801 lsolid5->SetVisAttributes(att("iron"));
802 G4Transform3D tsolid5_p1(G4RotateZ3D(M_PI / 16)*G4Translate3D(910, -45, Z0 - 20 - 15));
803 acs->AddPlacedVolume(lsolid5, tsolid5_p1);
804 // new G4PVPlacement(tsolid5_p1, lsolid5, "psolid5_p1", crystalSectorLogical, false, 0, overlap);
805 G4Transform3D tsolid5_p3(G4RotateZ3D(-M_PI / 16)*G4Translate3D(910, 45, Z0 - 20 - 15)*G4RotateZ3D(M_PI));
806 acs->AddPlacedVolume(lsolid5, tsolid5_p3);
807 // new G4PVPlacement(tsolid5_p3, lsolid5, "psolid5_p3", crystalSectorLogical, false, 0, overlap);
808
809 G4VSolid* solid7_p1 = new G4Box("fwd_solid7_p1", 130. / 2, 30. / 2, 40. / 2);
810 G4VSolid* solid7_p2 = new G4Box("fwd_solid7_p2", 131. / 2, 30. / 2, 40. / 2);
811 G4VSolid* solid7_p3 = new G4SubtractionSolid("fwd_solid7_p3", solid7_p1, solid7_p2, G4Transform3D(G4RotationMatrix(),
812 G4ThreeVector(0,
813 -3, 3)));
814 G4LogicalVolume* lsolid7 = new G4LogicalVolume(solid7_p3, Materials::get("SUS304"), "lsolid7", 0, 0, 0);
815 lsolid7->SetVisAttributes(att("iron"));
816 G4Transform3D tsolid7_p1(G4RotateZ3D(-M_PI / 16)*G4Translate3D(Ro - 8 - 65, 54, Z0 - 40. / 2)*G4RotateZ3D(M_PI));
817 acs->AddPlacedVolume(lsolid7, tsolid7_p1);
818 // new G4PVPlacement(tsolid7_p1, lsolid7, "psolid7_p1", crystalSectorLogical, false, 0, overlap);
819 G4Transform3D tsolid7_p2(G4RotateZ3D(M_PI / 16)*G4Translate3D(Ro - 8 - 65 - 8, -85, Z0 - 40. / 2));
820 acs->AddPlacedVolume(lsolid7, tsolid7_p2);
821 // new G4PVPlacement(tsolid7_p2, lsolid7, "psolid7_p2", crystalSectorLogical, false, 0, overlap);
822
823 // G4Transform3D tsolid7_p3(G4RotateZ3D(M_PI/16)*G4TranslateZ3D(1960 + 438 - 40./2 - 1)*G4TranslateX3D(1315)*G4TranslateY3D(-85));
824 // new G4PVPlacement(tsolid7_p3, lsolid7, "psolid7_p3", crystalSectorLogical, false, 0, overlap);
825 // G4Transform3D tsolid7_p4(G4RotateZ3D(M_PI/16)*G4TranslateZ3D(1960 + 438 - 40./2 - 1)*G4TranslateX3D(1315)*G4TranslateY3D(85)*G4RotateZ3D(M_PI));
826 // new G4PVPlacement(tsolid7_p4, lsolid7, "psolid7_p4", crystalSectorLogical, false, 0, overlap);
827
828 G4VSolid* solid8_p1 = new G4Box("fwd_solid8_p1", 120. / 2, 10. / 2, 42. / 2);
829 G4LogicalVolume* lsolid8 = new G4LogicalVolume(solid8_p1, Materials::get("SUS304"), "lsolid8", 0, 0, 0);
830 lsolid8->SetVisAttributes(att("iron"));
831 G4Transform3D tsolid8_p1(G4RotateZ3D(-M_PI / 16)*G4Translate3D(Ro - 8 - 60, 34, Z0 - 42. / 2));
832 acs->AddPlacedVolume(lsolid8, tsolid8_p1);
833 // new G4PVPlacement(tsolid8_p1, lsolid8, "psolid8_p1", crystalSectorLogical, false, 0, overlap);
834 // G4Transform3D tsolid8_p2(G4RotateZ3D(-M_PI/16)*G4TranslateZ3D(1960 + 438 - 42./2 - 1)*G4TranslateX3D(1320)*G4TranslateY3D(34)*G4RotateZ3D(M_PI));
835 // new G4PVPlacement(tsolid8_p2, lsolid8, "psolid8_p2", crystalSectorLogical, false, 0, overlap);
836
837
838 G4VSolid* solid9_p1 = new G4Box("fwd_solid9_p1", 26. / 2, 12. / 2, 26. / 2);
839 G4LogicalVolume* lsolid9 = new G4LogicalVolume(solid9_p1, Materials::get("SUS304"), "lsolid9", 0, 0, 0);
840 lsolid9->SetVisAttributes(att("iron"));
841 G4Transform3D tsolid9_p1(G4RotateZ3D(-M_PI / 16)*G4Translate3D(Ro - 8 - 120 + 26. / 2 + 17, 20 + 9 + 10 + 9 + 3,
842 Z0 - 37 + 26. / 2));
843 acs->AddPlacedVolume(lsolid9, tsolid9_p1);
844 // new G4PVPlacement(tsolid9_p1, lsolid9, "psolid9", crystalSectorLogical, false, 0, overlap);
845
846
847 G4VSolid* solid10_p1 = new G4Box("fwd_solid10_p1", 26. / 2, 12. / 2, 42. / 2);
848 G4LogicalVolume* lsolid10 = new G4LogicalVolume(solid10_p1, Materials::get("SUS304"), "lsolid10", 0, 0, 0);
849 lsolid10->SetVisAttributes(att("iron"));
850 // G4Transform3D tsolid5_p3(G4RotateZ3D(-M_PI/16)*G4Translate3D(910, 45, Z0-20-15)*G4RotateZ3D(M_PI));
851 G4Transform3D tsolid10_p1(G4RotateZ3D(-M_PI / 16)*G4Translate3D(910 + 650. / 2 - 221, 20 + 10 + 6 + 6, Z0 - 52 + 42. / 2));
852 acs->AddPlacedVolume(lsolid10, tsolid10_p1);
853 // new G4PVPlacement(tsolid10_p1, lsolid10, "psolid10", crystalSectorLogical, false, 0, overlap);
854
855 if (b_connectors) {
856 double t = 2, h20 = 32;
857 G4VSolid* solid_connector = new G4Box("fwd_solid_connector", (110 + 2 * 20) / 2, (250 + 2 * 20) / 2, h20 / 2);
858 G4VSolid* solid_connector2 = new G4Box("fwd_solid_connector2", (20 + t) / 2, (20 + t) / 2, 1.01 * h20 / 2);
859 G4VSolid* solid_connector3 = new G4Box("fwd_solid_connector3", 20 / 2, (250 + 2 * 20 + 2) / 2, h20 / 2);
860 solid_connector = new G4SubtractionSolid("fwd_solid_connector", solid_connector, solid_connector2, G4Translate3D(55 + (20 + t) / 2,
861 -140, 0));
862 solid_connector = new G4SubtractionSolid("fwd_solid_connector", solid_connector, solid_connector3, G4Translate3D(-70, 0, -t - 1));
863 G4LogicalVolume* lsolid_connector = new G4LogicalVolume(solid_connector, Materials::get("G4_AIR"), "lsolid_connector", 0,
864 0, 0);
865 lsolid_connector->SetVisAttributes(att("air"));
866 G4Transform3D tsolid_connector(G4Translate3D(1360 - 60, 0, Z0 - h20 / 2));
867 acs->AddPlacedVolume(lsolid_connector, tsolid_connector);
868 // new G4PVPlacement(tsolid_connector, lsolid_connector, "psolid20", crystalSectorLogical, false, 0, overlap);
869
870 auto lvolume = [&](int part, double dx, double dy, double dz) {
871 ostringstream ost(""); ost << "solid20_p" << part;
872 G4VSolid* sv = new G4Box(ost.str().c_str(), dx / 2, dy / 2, dz / 2);
873 ost.str(""); ost << "lsolid20_p" << part;
874 return new G4LogicalVolume(sv, Materials::get("A5052"), ost.str().c_str(), 0, 0, 0);
875 };
876
877 auto place = [&](G4LogicalVolume * lv, const G4Translate3D & move, int n) {
878 lv->SetVisAttributes(att("alum"));
879 ostringstream ost(""); ost << "phys_" << lv->GetName();
880 new G4PVPlacement(move, lv, ost.str().c_str(), lsolid_connector, false, n, overlap);
881 };
882 G4LogicalVolume* lv1 = lvolume(1, 20, 250 + 2 * 20, t);
883 place(lv1, G4Translate3D(-55 - 10, 0, h20 / 2 - t / 2), 0);
884 G4LogicalVolume* lv1_2 = lvolume(1, 20, 250 + 20, t);
885 place(lv1_2, G4Translate3D(55 + 10, 10, h20 / 2 - t / 2), 1);
886
887 G4LogicalVolume* lv2 = lvolume(2, 110., 20., t);
888 place(lv2, G4Translate3D(0, 250 / 2 + 10, h20 / 2 - t / 2), 0);
889 place(lv2, G4Translate3D(0, -250 / 2 - 10, h20 / 2 - t / 2), 1);
890
891 G4LogicalVolume* lv3 = lvolume(3, 110., t, h20 - t);
892 place(lv3, G4Translate3D(0, 250 / 2 + t / 2, -t / 2), 0);
893 place(lv3, G4Translate3D(0, -250 / 2 - t / 2, -t / 2), 1);
894
895 G4LogicalVolume* lv4 = lvolume(4, t, 250 + 2 * t, h20 - t);
896 place(lv4, G4Translate3D(55 + t / 2, 0, -t / 2), 0);
897 place(lv4, G4Translate3D(-55 - t / 2, 0, -t / 2), 1);
898
899 G4LogicalVolume* lv5 = lvolume(5, 7, 250, t);
900 place(lv5, G4Translate3D(55 - 7. / 2, 0, -h20 / 2 + t / 2 + t), 0);
901 place(lv5, G4Translate3D(-55 + 7. / 2, 0, -h20 / 2 + t / 2 + t), 1);
902
903 G4LogicalVolume* lv6 = lvolume(6, 110 - 14, 7, t);
904 place(lv6, G4Translate3D(0, 250 / 2 - 7. / 2, -h20 / 2 + t / 2 + t), 0);
905 place(lv6, G4Translate3D(0, -250 / 2 + 7. / 2, -h20 / 2 + t / 2 + t), 1);
906
907 G4LogicalVolume* lv7 = lvolume(7, 110, 250, t);
908 place(lv7, G4Translate3D(0, 0, -h20 / 2 + t / 2), 0);
909
910 // G4LogicalVolume *lv8 = lvolume(8, 90, 10, 30);
911 G4VSolid* p8_1 = new G4Box("fwd_solid20_p8_1", 90. / 2, 10. / 2, 30. / 2);
912 G4VSolid* p8_2 = new G4Box("fwd_solid20_p8_2", 88. / 2, 8. / 2, 30. / 2);
913 G4VSolid* sp8 = new G4SubtractionSolid("fwd_solid20_p8", p8_1, p8_2, G4TranslateZ3D(-1));
914 G4LogicalVolume* lv8 = new G4LogicalVolume(sp8, Materials::get("A5052"), "lsolid20_p8", 0, 0, 0);
915 lv8->SetVisAttributes(att("alum2"));
916 for (int i = 0; i < 10; i++) place(lv8, G4Translate3D(0, 25 * (i - 4.5), -h20 / 2 + t + 30 / 2), i);
917 }
918
919 if (b_boards) {
920 double hbv = 30;
921 G4VSolid* solid_board = new G4Box("fwd_solid_board", (210) / 2, (110) / 2, hbv / 2);
922 G4LogicalVolume* lsolid_board = new G4LogicalVolume(solid_board, Materials::get("G4_AIR"), "lsolid_board", 0, 0, 0);
923 lsolid_board->SetVisAttributes(att("air"));
924 for (int i = 0; i < 1; i++) {
925 // G4Transform3D t0 = G4RotateZ3D(M_PI/16*(1-2*i))*G4Translate3D(598-8+210/2, 0, Z0-52+8+5+hbv/2);
926 G4Transform3D t0 = G4Translate3D(598 - 8 + 210 / 2, 0, Z0 - 52 + 8 + 5 + hbv / 2);
927 G4Transform3D t1 = t0 * G4Translate3D(210, 110 / 2 + 5, 0);
928 G4Transform3D t2 = t0 * G4Translate3D(210, -110 / 2 + 5, 0);
929 G4Transform3D t3 = t0 * G4Translate3D(2 * 210, 110 / 2 + 5, 0);
930 G4Transform3D t4 = t0 * G4Translate3D(2 * 210, -110 / 2 + 5, 0);
931 // new G4PVPlacement(t0, lsolid_board, "phys_solid_board0", crystalSectorLogical, false, 0, overlap);
932 // new G4PVPlacement(t1, lsolid_board, "phys_solid_board1", crystalSectorLogical, false, 1, overlap);
933 // new G4PVPlacement(t2, lsolid_board, "phys_solid_board2", crystalSectorLogical, false, 2, overlap);
934 // new G4PVPlacement(t3, lsolid_board, "phys_solid_board3", crystalSectorLogical, false, 3, overlap);
935 // new G4PVPlacement(t4, lsolid_board, "phys_solid_board4", crystalSectorLogical, false, 4, overlap);
936 acs->AddPlacedVolume(lsolid_board, t0);
937 acs->AddPlacedVolume(lsolid_board, t1);
938 acs->AddPlacedVolume(lsolid_board, t2);
939 acs->AddPlacedVolume(lsolid_board, t3);
940 acs->AddPlacedVolume(lsolid_board, t4);
941 }
942 G4Material* boxmaterial = Materials::get("G4_GLASS_PLATE");
943 auto lvolumeb = [&](int part, double dx, double dy, double dz) {
944 ostringstream ost(""); ost << "fwd_sboard_p" << part;
945 G4VSolid* sv = new G4Box(ost.str().c_str(), dx / 2, dy / 2, dz / 2);
946 ost.str(""); ost << "lboard_p" << part;
947 return new G4LogicalVolume(sv, boxmaterial, ost.str().c_str(), 0, 0, 0);
948 };
949
950 const G4VisAttributes* asolid20 = att("plate");
951
952 auto placeb = [&](G4LogicalVolume * lv, const G4Translate3D & move, int n) {
953 lv->SetVisAttributes(asolid20);
954 ostringstream ost(""); ost << "phys_" << lv->GetName();
955 new G4PVPlacement(move, lv, ost.str().c_str(), lsolid_board, false, n, overlap);
956 };
957
958 double hboard = 2;
959 G4LogicalVolume* lb1 = lvolumeb(1, 210, 110, hboard);
960 placeb(lb1, G4Translate3D(0, 0, -hbv / 2 + hboard / 2), 0);
961
962 double wcon = 20, hcon = 40, hc = hbv - hboard;
963 G4VSolid* sv_connector_bundle = new G4Box("fwd_sv_connector_bundle", 4 * wcon / 2, hcon / 2, hc / 2);
964 G4LogicalVolume* lv_connector_bundle = new G4LogicalVolume(sv_connector_bundle, Materials::get("G4_AIR"),
965 "lv_connector_bundle", 0, 0, 0);
966 lv_connector_bundle->SetVisAttributes(att("air"));
967 new G4PVPlacement(G4Translate3D(-210 / 2 + 10 + wcon * 2, -110 / 2 + hcon / 2, -hbv / 2 + hboard + hc / 2), lv_connector_bundle,
968 "pv_connector_bundle", lsolid_board, false, 0, overlap);
969 new G4PVPlacement(G4Translate3D(-210 / 2 + 10 + wcon * 2, 10 + hcon / 2, -hbv / 2 + hboard + hc / 2), lv_connector_bundle,
970 "pv_connector_bundle", lsolid_board, false, 0, overlap);
971 new G4PVPlacement(G4Translate3D(10 + wcon * 2, -110 / 2 + hcon / 2, -hbv / 2 + hboard + hc / 2), lv_connector_bundle,
972 "pv_connector_bundle", lsolid_board, false, 0, overlap);
973 new G4PVPlacement(G4Translate3D(10 + wcon * 2, 10 + hcon / 2, -hbv / 2 + hboard + hc / 2), lv_connector_bundle,
974 "pv_connector_bundle", lsolid_board, false, 0, overlap);
975
976 G4VSolid* sv_crystal_connector = new G4Box("fwd_sv_crystal_connector", wcon / 2, hcon / 2, hc / 2);
977 G4LogicalVolume* lv_crystal_connector = new G4LogicalVolume(sv_crystal_connector, Materials::get("G4_AIR"),
978 "lv_crystal_connector", 0, 0, 0);
979 lv_crystal_connector->SetVisAttributes(att("air"));
980
981 new G4PVPlacement(G4Translate3D(-1.5 * 20, 0, 0), lv_crystal_connector, "pv_crystal_connector", lv_connector_bundle, false, 0,
982 overlap);
983 new G4PVPlacement(G4Translate3D(-0.5 * 20, 0, 0), lv_crystal_connector, "pv_crystal_connector", lv_connector_bundle, false, 1,
984 overlap);
985 new G4PVPlacement(G4Translate3D(0.5 * 20, 0, 0), lv_crystal_connector, "pv_crystal_connector", lv_connector_bundle, false, 2,
986 overlap);
987 new G4PVPlacement(G4Translate3D(1.5 * 20, 0, 0), lv_crystal_connector, "pv_crystal_connector", lv_connector_bundle, false, 3,
988 overlap);
989
990 G4VSolid* sv_crystal_connector_p1 = new G4Box("fwd_sv_crystal_connector_p1", 8 / 2, 30 / 2, 20 / 2);
991 G4LogicalVolume* lv_crystal_connector_p1 = new G4LogicalVolume(sv_crystal_connector_p1, Materials::get("SUS304"),
992 "lv_crystal_connector_p1", 0, 0, 0);
993 lv_crystal_connector_p1->SetVisAttributes(att("connector"));
994
995 new G4PVPlacement(G4Translate3D(-5, 0, -hc / 2 + 20. / 2), lv_crystal_connector_p1, "pv_crystal_connector_p1", lv_crystal_connector,
996 false, 0, overlap);
997
998 G4VSolid* sv_capacitor = new G4Tubs("fwd_sv_capacitor", 0, 5, 5. / 2, 0, 2 * M_PI);
999 G4LogicalVolume* lv_capacitor = new G4LogicalVolume(sv_capacitor, Materials::get("SUS304"), "lv_capacitor", 0, 0, 0);
1000 lv_capacitor->SetVisAttributes(att("capacitor"));
1001
1002 new G4PVPlacement(G4Translate3D(5, -15, -hc / 2 + 5. / 2), lv_capacitor, "pv_capacitor", lv_crystal_connector, false, 0, overlap);
1003 new G4PVPlacement(G4Translate3D(5, -5, -hc / 2 + 5. / 2), lv_capacitor, "pv_capacitor", lv_crystal_connector, false, 1, overlap);
1004 new G4PVPlacement(G4Translate3D(5, 5, -hc / 2 + 5. / 2), lv_capacitor, "pv_capacitor", lv_crystal_connector, false, 2, overlap);
1005 new G4PVPlacement(G4Translate3D(5, 15, -hc / 2 + 5. / 2), lv_capacitor, "pv_capacitor", lv_crystal_connector, false, 3, overlap);
1006
1007 G4VSolid* sv_board_connector_p1 = new G4Box("fwd_sv_board_connector_p1", 80. / 2, 8. / 2, 20. / 2);
1008 G4LogicalVolume* lv_board_connector_p1 = new G4LogicalVolume(sv_board_connector_p1, Materials::get("SUS304"),
1009 "lv_board_connector_p1", 0, 0, 0);
1010 lv_board_connector_p1->SetVisAttributes(att("connector"));
1011
1012 new G4PVPlacement(G4Translate3D(-210 / 2 + 10 + 80 / 2, 0, -hbv / 2 + hboard + 20. / 2), lv_board_connector_p1,
1013 "pv_board_connector_p1", lsolid_board, false, 0, overlap);
1014 new G4PVPlacement(G4Translate3D(210 / 2 - 10 - 80 / 2, 0, -hbv / 2 + hboard + 20. / 2), lv_board_connector_p1,
1015 "pv_board_connector_p1", lsolid_board, false, 1, overlap);
1016 }
1017
1018 G4Transform3D tr = G4RotateZ3D(M_PI / 16) * G4ReflectY3D();
1019 acs->MakeImprint(innervolumesector_logical, tr, 0, overlap);
1020 tr = G4RotateZ3D(-M_PI / 16);
1021 acs->MakeImprint(innervolumesector_logical, tr, 1, overlap);
1022 delete acs;
1023 }// end of ECL-004K102.pdf page 15
1024
1025 if (b_cover) {
1026 G4VSolid* solid8_p1 = new G4Tubs("fwd_solid8_p1", RI + tand(13.12) * (434 + 1) - 20 / cosd(13.12), 1415, 1. / 2, -M_PI / 16,
1027 M_PI / 8);
1028 G4VSolid* solid8_p2 = new G4Box("fwd_solid8_p2", 130. / 2, 270. / 2, 2. / 2);
1029 G4VSolid* solid8_p3 = new G4Tubs("fwd_solid8_p3", 0, 16, 2, 0, 2 * M_PI);
1030 G4VSolid* solid8_p4 = new G4Box("fwd_solid8_p4", 130. / 2, (75. - 2 * 16.) / 2, 2. / 2);
1031 double width_p5 = 180;
1032 G4VSolid* solid8_p5 = new G4Box("fwd_solid8_p5", width_p5 / 2, 2.5 + 75. / 2, 2. / 2);
1033
1034 double xx0 = 1415 - 47.8715;
1035 G4VSolid* solid8 = new G4SubtractionSolid("fwd_solid8", solid8_p1, solid8_p2, G4TranslateX3D(xx0 - 130. / 2));
1036 double xx1 = xx0 + 1.7;
1037 solid8 = new G4SubtractionSolid("fwd_solid8", solid8, solid8_p3, G4Translate3D(xx1, 159.5, 0));
1038 solid8 = new G4SubtractionSolid("fwd_solid8", solid8, solid8_p3, G4Translate3D(xx1, 202.5, 0));
1039 solid8 = new G4SubtractionSolid("fwd_solid8", solid8, solid8_p4, G4Translate3D(xx1 - 16 + 130. / 2, (202.5 + 159.5) / 2, 0));
1040 solid8 = new G4SubtractionSolid("fwd_solid8", solid8, solid8_p5, G4Translate3D(xx1 + width_p5 / 2, (202.5 + 159.5) / 2, 0));
1041 solid8 = new G4SubtractionSolid("fwd_solid8", solid8, solid8_p5, G4Translate3D(xx1 - 130 + (1230.77 - 1230.88), -177.57,
1042 0)*G4RotateZ3D(-M_PI / 16)*G4Translate3D(width_p5 / 2, -75. / 2, 0));
1043
1044 // for(int i=0;i<100000;i++){
1045 // G4ThreeVector v = solid8->GetPointOnSurface();
1046 // G4cout<<v.x()<<" "<<v.y()<<" "<<v.z()<<"\n";
1047 // }
1048
1049 G4LogicalVolume* lsolid8 = new G4LogicalVolume(solid8, Materials::get("A5052"), "lsolid8", 0, 0, 0);
1050 lsolid8->SetVisAttributes(att("alum"));
1051 for (int i = 0; i < 8; i++) {
1052 G4Transform3D tc = gTrans * G4Translate3D(0, 0, 3 + 434 + 0.5) * G4RotateZ3D(M_PI / 8 + i * M_PI / 4);
1053 new G4PVPlacement(tc * G4RotateZ3D(M_PI / 16), lsolid8, suf("ECL_Forward_cover", 0 + 2 * i), top, false, 0 + 2 * i, overlap);
1054 G4ReflectionFactory::Instance()->Place(tc * G4RotateZ3D(-M_PI / 16)*G4ReflectY3D(), suf("ECL_Forward_cover", 0 + 2 * i), lsolid8,
1055 top, false,
1056 1 + 2 * i, overlap);
1057 }
1058 }
1059 // end of ECL-004K102.pdf page 11 - 12
1060
1061}

◆ get_pa_box_height()

double get_pa_box_height ( ) const
inlineprivate

Getter for preamplifier box height (hard-coded to be 2)

Definition at line 88 of file GeoECLCreator.h.

88{return 2;}

◆ get_preamp()

G4LogicalVolume * get_preamp ( ) const
private

Get Logical volume of preamplifier.

Definition at line 148 of file GeoECLCreator.cc.

149{
150 static G4LogicalVolume* lv_preamplifier = nullptr;
151 if (lv_preamplifier == nullptr) {
152 G4VSolid* sv_preamplifier = new G4Box("sv_preamplifier", 58. / 2, 51. / 2, get_pa_box_height() / 2);
153 lv_preamplifier = new G4LogicalVolume(sv_preamplifier, Materials::get("A5052"), "lv_preamplifier", 0, 0, 0);
154 G4VSolid* sv_diode = new G4Box("sv_diode", 20. / 2, 20. / 2, 0.3 / 2);
155 G4LogicalVolume* lv_diode = new G4LogicalVolume(sv_diode, Materials::get("G4_Si"), "lv_diode", 0, 0, 0);
156 lv_diode->SetUserLimits(new G4UserLimits(0.01));
157 lv_diode->SetSensitiveDetector(m_sensediode);
158 new G4PVPlacement(G4TranslateZ3D(-get_pa_box_height() / 2 + 0.3 / 2), lv_diode, "pv_diode", lv_preamplifier, false, 0, m_overlap);
159
160 // G4VSolid* sv_diode = new G4Box("sv_diode", 10. / 2, 20. / 2, 0.3 / 2);
161 // G4LogicalVolume* lv_diode = new G4LogicalVolume(sv_diode, Materials::get("G4_Si"), "lv_diode", 0, 0, 0);
162 // lv_diode->SetUserLimits(new G4UserLimits(0.01));
163 // new G4PVPlacement(G4Translate3D(-5, 0, -pa_box_height / 2 + 0.3 / 2), lv_diode, "pv_diode1", lv_preamplifier, false, 1, overlap);
164 // new G4PVPlacement(G4Translate3D(5, 0, -pa_box_height / 2 + 0.3 / 2), lv_diode, "pv_diode2", lv_preamplifier, false, 2, overlap);
165
166 lv_preamplifier->SetVisAttributes(att("preamp"));
167 }
168 return lv_preamplifier;
169}

◆ wrapped_crystal()

G4LogicalVolume * wrapped_crystal ( const shape_t s,
const std::string &  endcap,
double  wrapthickness 
)
private

Wrapped crystal.

Definition at line 90 of file GeoECLCreator.cc.

91{
92 std::string prefix("sv_"); prefix += endcap; prefix += "_wrap";
93 G4Translate3D tw;
94 G4VSolid* wrapped_crystal = s->get_solid(prefix, wrapthickness, tw);
95 std::string name("lv_"); name += endcap + "_wrap_" + std::to_string(s->nshape);
96 G4Material* wrap = nullptr;
97 if (wrapthickness < 0.170)
98 wrap = Materials::get("WRAP170");
99 else if (wrapthickness < 0.200)
100 wrap = Materials::get("WRAP200");
101 else
102 wrap = Materials::get("WRAP250");
103 G4LogicalVolume* wrapped_logical = new G4LogicalVolume(wrapped_crystal, wrap, name.c_str(), 0, 0, 0);
104 wrapped_logical->SetVisAttributes(att("wrap"));
105
106 prefix = "sv_"; prefix += endcap; prefix += "_crystal";
107 G4Translate3D tc;
108 G4VSolid* crystal_solid = s->get_solid(prefix, 0, tc);
109 name = "lv_" + endcap + "_crystal_" + std::to_string(s->nshape);
110 G4LogicalVolume* crystal_logical = new G4LogicalVolume(crystal_solid, Materials::get("G4_CESIUM_IODIDE"), name.c_str(),
111 0, 0, 0);
112 crystal_logical->SetVisAttributes(att("cryst"));
113 crystal_logical->SetSensitiveDetector(m_sensitive);
114
115 new G4PVPlacement(nullptr, G4ThreeVector(), crystal_logical, name.c_str(), wrapped_logical, false, 0, 0);
116 return wrapped_logical;
117}

Member Data Documentation

◆ m_atts

std::map<std::string, G4VisAttributes*> m_atts
private

Vector of background-Sensitive detectors.

Definition at line 98 of file GeoECLCreator.h.

◆ m_overlap

int m_overlap
private

overlap

Definition at line 100 of file GeoECLCreator.h.

◆ m_sap

const ECLCrystalsShapeAndPosition* m_sap
private

pointer to a storage with crystal shapes and positions

Definition at line 91 of file GeoECLCreator.h.

◆ m_sensediode

Simulation::SensitiveDetectorBase* m_sensediode
private

Sensitive diode.

Definition at line 96 of file GeoECLCreator.h.

◆ m_sensitive

Simulation::SensitiveDetectorBase* m_sensitive
private

Sensitive detector.

Definition at line 94 of file GeoECLCreator.h.


The documentation for this class was generated from the following files: