Belle II Software  release-08-01-10
GeoARICHCreator.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <arich/geometry/GeoARICHCreator.h>
10 #include <arich/dbobjects/tessellatedSolidStr.h>
11 
12 #include <geometry/Materials.h>
13 #include <geometry/CreatorFactory.h>
14 #include <geometry/utilities.h>
15 #include <framework/gearbox/Unit.h>
16 #include <framework/logging/Logger.h>
17 #include <arich/simulation/SensitiveDetector.h>
18 #include <arich/simulation/SensitiveAero.h>
19 #include <simulation/background/BkgSensitiveDetector.h>
20 #include <arich/dbobjects/ARICHPositionElement.h>
21 
22 #include <cmath>
23 #include <boost/format.hpp>
24 #include <boost/foreach.hpp>
25 #include <boost/algorithm/string.hpp>
26 
27 // Geant4
28 #include <G4LogicalVolume.hh>
29 #include <G4PVPlacement.hh>
30 #include <G4AssemblyVolume.hh>
31 #include <G4UnionSolid.hh>
32 #include <G4LogicalSkinSurface.hh>
33 #include <G4OpticalSurface.hh>
34 
35 // Geant4 Shapes
36 #include <G4Box.hh>
37 #include <G4Tubs.hh>
38 #include <G4Trap.hh>
39 #include <G4Torus.hh>
40 #include <G4TessellatedSolid.hh>
41 #include <G4TriangularFacet.hh>
42 
43 #include <G4SubtractionSolid.hh>
44 #include <G4Region.hh>
45 #include <G4Material.hh>
46 #include <TVector3.h>
47 
48 using namespace std;
49 using namespace boost;
50 using namespace CLHEP;
51 
52 namespace Belle2 {
58  using namespace geometry;
59 
60  namespace arich {
61 
62  //-----------------------------------------------------------------
63  // Register the Creator
64  //-----------------------------------------------------------------
65 
66  geometry::CreatorFactory<GeoARICHCreator> GeoARICHFactory("ARICHCreator");
67 
68  //-----------------------------------------------------------------
69  // Implementation
70  //-----------------------------------------------------------------
71 
72  GeoARICHCreator::GeoARICHCreator(): m_isBeamBkgStudy(0)
73  {
74  m_sensitive = NULL;
75  m_sensitiveAero = NULL;
76  }
77 
79  {
80  delete m_sensitive;
81  delete m_sensitiveAero;
82  G4LogicalSkinSurface::CleanSurfaceTable();
83 
84  }
85 
86 
87  void GeoARICHCreator::createGeometry(G4LogicalVolume& topVolume, GeometryTypes)
88  {
89 
93 
94  // print geometry configuration
95  // m_config.print();
96 
97  //Build envelope
98  G4Tubs* envelopeTube = new G4Tubs("Envelope", m_config.getMasterVolume().getInnerRadius(),
100  G4Material* material = Materials::get(m_config.getMasterVolume().getMaterial());
101 
102  // check of material and its refractive index
103  if (!material) { B2FATAL("Material ARICH_Air required for ARICH master volume could not be found");}
104  if (!getAvgRINDEX(material))
105  B2WARNING("Material ARICH_Air required by ARICH master volume has no specified refractive index. Continuing, but no photons in ARICH will be propagated.");
106 
107  // create and place
108  G4LogicalVolume* masterLV = new G4LogicalVolume(envelopeTube, material, "ARICH.masterVolume");
109  setVisibility(*masterLV, false);
110 
111  // Set up region for production cuts
112  G4Region* aRegion = new G4Region("ARICHEnvelope");
113  masterLV->SetRegion(aRegion);
114  aRegion->AddRootLogicalVolume(masterLV);
115 
116  G4RotationMatrix rotMaster;
117  double rot_x = m_config.getMasterVolume().getRotationX();
118  double rot_y = m_config.getMasterVolume().getRotationY();
119  double rot_z = m_config.getMasterVolume().getRotationZ();
120  double tr_x = m_config.getMasterVolume().getPosition().X();
121  double tr_y = m_config.getMasterVolume().getPosition().Y();
122  double tr_z = m_config.getMasterVolume().getPosition().Z();
123 
124  if (m_config.useGlobalDisplacement()) {
128  tr_x += m_config.getGlobalDisplacement().getX();
129  tr_y += m_config.getGlobalDisplacement().getY();
130  tr_z += m_config.getGlobalDisplacement().getZ();
131  B2WARNING("ARICH global displacement parameters from DB will be taken into account.");
132  }
133  rotMaster.rotateX(rot_x);
134  rotMaster.rotateY(rot_y);
135  rotMaster.rotateZ(rot_z);
136 
137  G4ThreeVector transMaster(tr_x, tr_y, tr_z);
138 
139  new G4PVPlacement(G4Transform3D(rotMaster, transMaster), masterLV, "ARICH.MasterVolume", &topVolume, false, 1);
140 
142 
143  // Z shift for placing volumes in ARICH master volume
144  double zShift = 0;
145 
146  // build photon detector logical volume
147  G4LogicalVolume* detPlaneLV = buildDetectorPlane(m_config);
148  G4LogicalVolume* detSupportPlateLV = buildDetectorSupportPlate(m_config);
149 
150  // Build merger PCB logical volume
151  G4LogicalVolume* mergerLV = buildMergerPCBEnvelopePlane(m_config);
152 
153  // Build the cables envelop with effective material describing cables
154  G4LogicalVolume* cablesLV = buildCables(m_config.getCablesEnvelope());
155 
156  // Build cooling system assembly envelope plane
157  //G4LogicalVolume* coolingLV = buildCoolingEnvelopePlane(m_config.getCoolingGeometry());
158 
159  // Build aerogel logical volume
160  G4LogicalVolume* aeroPlaneLV;
161  if (!m_config.getAerogelPlane().isSimple()) {
162  aeroPlaneLV = buildAerogelPlane(m_config);
163  } else {
164  B2INFO("GeoARICHCreator: using simple cosmic test geometry.");
165  aeroPlaneLV = buildSimpleAerogelPlane(m_config);
166  }
167 
168  G4RotationMatrix rotDetPlane;
169  rotDetPlane.rotateX(m_config.getDetectorPlane().getRotationX());
170  rotDetPlane.rotateY(m_config.getDetectorPlane().getRotationY());
171  rotDetPlane.rotateZ(m_config.getDetectorPlane().getRotationZ());
172 
173  //
174  // For the moment rotation of merger PCB board envelope is not foreseen.
175  //
176  G4RotationMatrix rotMergerPlane;
177  rotMergerPlane.rotateX(0.0);
178  rotMergerPlane.rotateY(0.0);
179  rotMergerPlane.rotateZ(0.0);
180 
181  //
182  // For the moment rotation of cables envelope is not foreseen.
183  //
184  G4RotationMatrix rotCablesPlane;
185  rotCablesPlane.rotateX(0.0);
186  rotCablesPlane.rotateY(0.0);
187  rotCablesPlane.rotateZ(0.0);
188 
189  G4RotationMatrix rotCoolingPlane;
190  rotCoolingPlane.rotateX(0.0);
191  rotCoolingPlane.rotateY(0.0);
192  rotCoolingPlane.rotateZ(0.0);
193 
194  G4RotationMatrix rotAeroPlane;
195  rotAeroPlane.rotateX(m_config.getAerogelPlane().getRotationX());
196  rotAeroPlane.rotateY(m_config.getAerogelPlane().getRotationY());
197  rotAeroPlane.rotateZ(m_config.getAerogelPlane().getRotationZ());
198 
199  G4ThreeVector transDetPlane(m_config.getDetectorPlane().getPosition().X(), m_config.getDetectorPlane().getPosition().Y(),
200  m_config.getDetectorPlane().getPosition().Z() + zShift);
201 
202  G4ThreeVector transDetSupportPlate(m_config.getDetectorPlane().getPosition().X(), m_config.getDetectorPlane().getPosition().Y(),
204 
205  G4ThreeVector transMergerPlate(m_config.getMergerGeometry().getEnvelopeCenterPosition().X() * mm,
208 
209  G4ThreeVector transCablesPlate(m_config.getCablesEnvelope().getEnvelopeCenterPosition().X(),
212 
213  G4ThreeVector transCoolingPlate(m_config.getCoolingGeometry().getEnvelopeCenterPosition().X() * mm,
216 
217  G4ThreeVector transAeroPlane(m_config.getAerogelPlane().getPosition().X(),
219  m_config.getAerogelPlane().getPosition().Z() + zShift);
220 
221 
222  new G4PVPlacement(G4Transform3D(rotDetPlane, transDetPlane), detPlaneLV, "ARICH.detPlane", masterLV, false, 1);
223  new G4PVPlacement(G4Transform3D(rotDetPlane, transDetSupportPlate), detSupportPlateLV, "ARICH.detSupportPlane", masterLV, false, 1);
224  if (mergerLV) new G4PVPlacement(G4Transform3D(rotMergerPlane, transMergerPlate), mergerLV, "ARICH.mergerPlane", masterLV, false, 1);
225  new G4PVPlacement(G4Transform3D(rotCablesPlane, transCablesPlate), cablesLV, "ARICH.cablesPlane", masterLV, false, 1);
226  //new G4PVPlacement(G4Transform3D(rotCoolingPlane, transCoolingPlate), coolingLV, "ARICH.coolingPlane", masterLV, false, 1);
227  new G4PVPlacement(G4Transform3D(rotAeroPlane, transAeroPlane), aeroPlaneLV, "ARICH.aeroPlane", masterLV, false, 1);
228 
229  /*
230  // build and place cooling test plates
231  G4LogicalVolume* coolingTestPlatesLV = buildCoolingTestPlate(m_config.getCoolingGeometry());
232  double coolingTestPlateAngle = 0.0;
233  double coolingTestPlateR = 0.0;
234  for (unsigned i = 0; i < m_config.getCoolingGeometry().getCoolingTestPlatePosR().size(); i++) {
235  coolingTestPlateAngle = m_config.getCoolingGeometry().getCoolingTestPlatePosPhi().at(i) * deg;
236  coolingTestPlateR = m_config.getCoolingGeometry().getCoolingTestPlatePosR().at(i) * mm;
237  G4RotationMatrix rotCoolingTestPlate;
238  rotCoolingTestPlate.rotateZ(coolingTestPlateAngle);
239  G4ThreeVector transCoolingTestPlate(coolingTestPlateR * cos(coolingTestPlateAngle),
240  coolingTestPlateR * sin(coolingTestPlateAngle),
241  m_config.getCoolingGeometry().getCoolingTestPlatePosZ0().at(i) * mm + zShift);
242  new G4PVPlacement(G4Transform3D(rotCoolingTestPlate, transCoolingTestPlate), coolingTestPlatesLV, "ARICH.coolingTestPlates",
243  masterLV, false, i);
244  }
245  */
246 
247  // build and place mirrors
248  G4LogicalVolume* mirrorLV = buildMirror(m_config);
249  int nMirrors = m_config.getMirrors().getNMirrors();
250  int mirStart = m_config.getMirrors().getStartAngle();
251  int mirZPos = m_config.getMirrors().getZPosition();
252  int mirRad = m_config.getMirrors().getRadius();
253 
254  double angl = mirStart;
255  double dphi = 2 * M_PI / nMirrors;
256  if (m_config.useMirrorDisplacement()) B2WARNING("ARICH mirrors displacement parameters from DB will be used.");
257  for (int i = 1; i < nMirrors + 1; i++) {
258  double mrot_x = 0;
259  double mrot_y = 0;
260  double mrot_z = angl;
261  double mtr_x = mirRad * cos(angl);
262  double mtr_y = mirRad * sin(angl);
263  double mtr_z = mirZPos + zShift;
264  if (m_config.useMirrorDisplacement()) {
266  mrot_x += displ.getAlpha();
267  mrot_y += displ.getBeta();
268  mrot_z += displ.getGamma();
269  mtr_x += displ.getX();
270  mtr_y += displ.getY();
271  mtr_z += displ.getZ();
272  }
273  G4RotationMatrix rotMirror;
274  rotMirror.rotateX(mrot_x);
275  rotMirror.rotateY(mrot_y);
276  rotMirror.rotateZ(mrot_z);
277  G4ThreeVector transMirror(mtr_x, mtr_y, mtr_z);
278  new G4PVPlacement(G4Transform3D(rotMirror, transMirror), mirrorLV, "ARICH.mirrorPlate", masterLV, false, i);
279  angl += dphi;
280  }
281 
282  // build additional components of support structure
283 
284  const ARICHGeoSupport& supportPar = m_config.getSupportStructure();
285  std::vector<G4LogicalVolume*> shieldLV(2);
286  std::vector<double> shieldZ(2);
287  int nTubes = supportPar.getNTubes();
288  for (int i = 0; i < nTubes; i++) {
289  G4Tubs* tube = new G4Tubs(supportPar.getTubeName(i), supportPar.getTubeInnerR(i), supportPar.getTubeOuterR(i),
290  supportPar.getTubeLength(i) / 2., 0, 2.*M_PI);
291  G4Material* tubeMaterial = Materials::get(supportPar.getTubeMaterial(i));
292 
293  if (supportPar.getTubeName(i) == "ARICH.NeutronShield1") {
294  shieldLV[0] = new G4LogicalVolume(tube, tubeMaterial, supportPar.getTubeName(i));
295  shieldZ[0] = supportPar.getTubeZPosition(i) + supportPar.getTubeLength(i) / 2.;
296  continue;
297  } else if (supportPar.getTubeName(i) == "ARICH.NeutronShield2") {
298  shieldLV[1] = new G4LogicalVolume(tube, tubeMaterial, supportPar.getTubeName(i));
299  shieldZ[1] = supportPar.getTubeZPosition(i) + supportPar.getTubeLength(i) / 2.;
300  continue;
301  }
302 
303  G4LogicalVolume* tubeLV = new G4LogicalVolume(tube, tubeMaterial, supportPar.getTubeName(i));
304 
305  new G4PVPlacement(G4Transform3D(G4RotationMatrix(), G4ThreeVector(0, 0,
306  supportPar.getTubeZPosition(i) + supportPar.getTubeLength(i) / 2. + zShift)), tubeLV, supportPar.getTubeName(i), masterLV, false,
307  1);
308  }
309 
310  const std::vector<double> wedge1 = supportPar.getWedge(1);
311  const std::vector<double> wedge2 = supportPar.getWedge(2);
312  G4Material* wedge1Material = Materials::get(supportPar.getWedgeMaterial(1));
313  G4Material* wedge2Material = Materials::get(supportPar.getWedgeMaterial(2));
314  G4AssemblyVolume* assemblyWedge1 = new G4AssemblyVolume();
315  G4AssemblyVolume* assemblyWedge2 = new G4AssemblyVolume();
316  makeJoint(wedge1Material, wedge1, assemblyWedge1);
317  makeJoint(wedge2Material, wedge2, assemblyWedge2);
318 
319  G4Transform3D Tr;
320  int nWedge = supportPar.getNWedges();
321  for (int i = 0; i < nWedge; i++) {
322  G4RotationMatrix Rx;
323  int wgtype = supportPar.getWedgeType(i);
324  double phi = supportPar.getWedgePhi(i);
325  Rx.rotateZ(phi);
326  double r = supportPar.getWedgeR(i);
327  double z = supportPar.getWedgeZ(i);
328  G4ThreeVector Tb(r * cos(phi), r * sin(phi), z);
329 
330  if (shieldLV[0]) {
331  z -= shieldZ[0];
332  Tb.setZ(z);
333  Tr = G4Transform3D(Rx, Tb);
334  if (wgtype == 1) assemblyWedge1->MakeImprint(shieldLV[0], Tr);
335  else if (wgtype == 2) assemblyWedge2->MakeImprint(shieldLV[0], Tr);
336  continue;
337  }
338 
339  Tr = G4Transform3D(Rx, Tb);
340  if (wgtype == 1) assemblyWedge1->MakeImprint(masterLV, Tr);
341  else if (wgtype == 2) assemblyWedge2->MakeImprint(masterLV, Tr);
342  else B2ERROR("GeoARICHCreator: invalid support wedge type!");
343  }
344 
345  // place neutron shield volumes
346  if (shieldLV[0])
347  new G4PVPlacement(G4Transform3D(G4RotationMatrix(), G4ThreeVector(0, 0, shieldZ[0])), shieldLV[0], "ARICH.NeutronShield1", masterLV,
348  false, 1);
349  if (shieldLV[1])
350  new G4PVPlacement(G4Transform3D(G4RotationMatrix(), G4ThreeVector(0, 0, shieldZ[1])), shieldLV[1], "ARICH.NeutronShield2", masterLV,
351  false, 1);
352 
353  // place mirror holders, for now skip
354  /*
355  double mirSupX = 17.;
356  double mirSupY = 14.;
357  double mirSupLen = 131;
358  double mirSupLen1 = 201;
359 
360  std::vector<G4TwoVector> polygon;
361  polygon.assign(12, G4TwoVector());
362  polygon[11] = G4TwoVector(-mirSupX / 2., -mirSupY / 2.);
363  polygon[10] = G4TwoVector(-mirSupX / 2., 2.);
364  polygon[9] = G4TwoVector(-mirSupX / 4., 2. - 0.75);
365  polygon[8] = G4TwoVector(-mirSupX / 4., 5. - 0.75);
366  polygon[7] = G4TwoVector(-mirSupX / 2., 5);
367  polygon[6] = G4TwoVector(-mirSupX / 2., mirSupY / 2.);
368  polygon[5] = G4TwoVector(mirSupX / 2., mirSupY / 2.);
369  polygon[4] = G4TwoVector(mirSupX / 2., 6);
370  polygon[3] = G4TwoVector(mirSupX / 4., 6. - 0.75);
371  polygon[2] = G4TwoVector(mirSupX / 4., 2. - 0.75);
372  polygon[1] = G4TwoVector(mirSupX / 2., 2.);
373  polygon[0] = G4TwoVector(mirSupX / 2., -mirSupY / 2.);
374 
375  std::vector<G4ExtrudedSolid::ZSection> zsections;
376  zsections.push_back(G4ExtrudedSolid::ZSection(-mirSupLen / 2., G4TwoVector(0, 0), 1));
377  zsections.push_back(G4ExtrudedSolid::ZSection(mirSupLen / 2., G4TwoVector(0, 0), 1));
378 
379  G4ExtrudedSolid* shape = new G4ExtrudedSolid("bla", polygon, zsections);
380  G4LogicalVolume* shapeLV = new G4LogicalVolume(shape, wedge1Material, "mirrorsupport");
381 
382  G4Box* shape1 = new G4Box("shape1", mirSupX / 2., 2. / 2., mirSupLen1 / 2.);
383  G4LogicalVolume* shape1LV = new G4LogicalVolume(shape1, wedge1Material, "mirrorholder");
384 
385  G4AssemblyVolume* mirroHolder = new G4AssemblyVolume();
386  G4RotationMatrix Rh;
387  G4ThreeVector Th(0, 0, 0);
388  Tr = G4Transform3D(Rh, Th);
389  mirroHolder->AddPlacedVolume(shapeLV, Tr);
390 
391  Th.setZ(-35.0);
392  Th.setY(-mirSupY / 2. - 1.);
393  Tr = G4Transform3D(Rh, Th);
394  mirroHolder->AddPlacedVolume(shape1LV, Tr);
395 
396  double rholder = m_config.getDetectorPlane().getSupportOuterR() - mirSupY / 2.;
397 
398  angl = dphi / 2.;
399  for (int i = 1; i < nMirrors + 1; i++) {
400  G4RotationMatrix rotMirror;
401  rotMirror.rotateX(M_PI);
402  rotMirror.rotateZ(-M_PI / 2. + angl);
403  G4ThreeVector transMirror(rholder * cos(angl), rholder * sin(angl), mirZPos + zShift);
404  Tr = G4Transform3D(rotMirror, transMirror);
405  mirroHolder->MakeImprint(masterLV, Tr);
406  angl += dphi;
407  }
408  */
409 
410  // temporary for simple cosmic test
411  // if using simple configuration, place scinitilators
413  int nBoxes = supportPar.getNBoxes();
414  for (int i = 0; i < nBoxes; i++) {
415  ARICHGeoSupport::box box = supportPar.getBox(i);
416  G4Box* scintBox = new G4Box("scintBox", box.size[0] * 10. / 2., box.size[1] * 10. / 2., box.size[2] * 10. / 2.);
417  G4Material* scintMaterial = Materials::get(box.material);
418  G4LogicalVolume* scintLV = new G4LogicalVolume(scintBox, scintMaterial, box.name);
419  scintLV->SetSensitiveDetector(m_sensitiveAero);
420  G4RotationMatrix rotScint;
421  rotScint.rotateX(box.rotation[0]);
422  rotScint.rotateY(box.rotation[1]);
423  rotScint.rotateZ(box.rotation[2]);
424  TVector3 transScintTV(box.position[0], box.position[1], box.position[2]);
425  transScintTV = m_config.getMasterVolume().pointToGlobal(transScintTV);
426  B2INFO("GeoARICHCreator: Scintilator " << box.name << " placed at global: " << transScintTV.X() << " " << transScintTV.Y() << " " <<
427  transScintTV.Z());
428  G4ThreeVector transScint(box.position[0] * 10., box.position[1] * 10., box.position[2] * 10.);
429  new G4PVPlacement(G4Transform3D(rotScint, transScint), scintLV, "scintilator", masterLV, false, 1);
430  }
431  }
432 
434 
435  delete assemblyWedge1;
436  delete assemblyWedge2;
437 
438  return;
439 
440  }
441 
442 
444  {
445 
446  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
447 
448  const std::vector<double>& params = aeroGeo.getSimpleParams();
449 
450  // support plane
451  double rin = aeroGeo.getSupportInnerR();
452  double rout = aeroGeo.getSupportOuterR();
453  double thick = aeroGeo.getSupportThickness();
454  string supportMat = aeroGeo.getSupportMaterial();
455  double wallHeight = aeroGeo.getWallHeight();
456  G4Material* supportMaterial = Materials::get(supportMat);
457  G4Material* gapMaterial = Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
458 
459  // master volume
460  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight) / 2., 0, 2 * M_PI);
461  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
462 
463  // support plate
464  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
465  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
466  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
467 
468  unsigned nLayer = aeroGeo.getNLayers();
469 
470  // loop over layers
471  double zLayer = 0;
472  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
473  double layerThick = aeroGeo.getLayerThickness(iLayer);
474  std::stringstream tileName;
475  tileName << "aerogelTile_" << iLayer;
476 
477  G4Box* tileShape = new G4Box(tileName.str(), params[0] * 10. / 2., params[1] * 10. / 2., layerThick / 2.);
478 
479  G4Material* aeroMaterial = Materials::get(aeroGeo.getLayerMaterial(iLayer));
480 
481  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, aeroMaterial, string("ARICH.") + tileName.str());
482 
483  G4ThreeVector transTile(params[2] * 10., params[3] * 10., (thick + layerThick - wallHeight) / 2. + zLayer);
484  G4RotationMatrix Ra;
485  Ra.rotateZ(params[4]);
486  new G4PVPlacement(G4Transform3D(Ra, transTile), tileLV, string("ARICH.") + tileName.str(), aerogelPlaneLV, false, iLayer);
487 
488  zLayer += layerThick;
489  }
490 
491  new G4PVPlacement(G4Translate3D(0., 0., -wallHeight / 2.), supportTubeLV, "ARICH.AerogelSupportPlate", aerogelPlaneLV, false, 1);
492 
493  return aerogelPlaneLV;
494  }
495 
496  G4LogicalVolume* GeoARICHCreator::buildAerogelPlane(const ARICHGeometryConfig& detectorGeo)
497  {
498  if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 0)
499  return buildAerogelPlaneAveragedOverLayers(detectorGeo);
500  else if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 1)
501  return buildAerogelPlaneWithIndividualTilesProp(detectorGeo);
502  else
503  B2ERROR("GeoARICHCreator::buildAerogelPlane --> getFullAerogelMaterialDescriptionKey() is wrong");
504  return NULL;
505  }
506 
508  {
509 
510  //cout<<"GeoARICHCreator::buildAerogelPlaneAveragedOverLayers(const ARICHGeometryConfig& detectorGeo)"<<endl;
511 
512  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
513 
514  // support plane
515  double rin = aeroGeo.getSupportInnerR();
516  double rout = aeroGeo.getSupportOuterR();
517  double thick = aeroGeo.getSupportThickness();
518  double wallThick = aeroGeo.getWallThickness();
519  double wallHeight = aeroGeo.getWallHeight();
520  string supportMat = aeroGeo.getSupportMaterial();
521  G4Material* supportMaterial = Materials::get(supportMat);
522  G4Material* gapMaterial = Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
523  G4Material* imgMaterial = Materials::get("ARICH_Air");
524  // master volume
525 
526  //double imgTubeLen = 0.5; // if changed, change position of aerogel plane also in main function!!
527  double imgTubeLen = aeroGeo.getImgTubeThickness();
528  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight + imgTubeLen) / 2., 0, 2 * M_PI);
529  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
530 
531  // support plate
532  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
533  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
534  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
535 
536  // imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
537  G4Tubs* imgTube = new G4Tubs("imgTube", rin, rout, imgTubeLen / 2., 0, 2 * M_PI);
538  G4LogicalVolume* imgTubeLV = new G4LogicalVolume(imgTube, imgMaterial, "ARICH.AerogelImgPlate");
539  imgTubeLV->SetSensitiveDetector(m_sensitiveAero);
540 
541 
542  // read radiuses of aerogel slot aluminum walls
543  std::vector<double> wallR;
544  unsigned nRing = aeroGeo.getNRings();
545  for (unsigned iRing = 1; iRing < nRing + 1; iRing++) {
546  wallR.push_back(aeroGeo.getRingRadius(iRing));
547  }
548 
549  unsigned nLayer = aeroGeo.getNLayers();
550  double tileGap = aeroGeo.getTileGap();
551  G4Transform3D transform = G4Translate3D(0., 0., (thick - imgTubeLen) / 2.);
552 
553  for (unsigned iRing = 0; iRing < nRing; iRing++) {
554 
555  // aluminum walls between tile rings (r wall)
556  std::stringstream wallName;
557  wallName << "supportWallR_" << iRing + 1;
558  G4Tubs* supportWall = new G4Tubs(wallName.str().c_str(), wallR[iRing], wallR[iRing] + wallThick, wallHeight / 2., 0, 2 * M_PI);
559  G4LogicalVolume* supportWallLV = new G4LogicalVolume(supportWall, supportMaterial, string("ARICH.") + wallName.str().c_str());
560 
561  new G4PVPlacement(transform, supportWallLV, string("ARICH.") + wallName.str().c_str(), aerogelPlaneLV, false, 0);
562 
563  if (iRing == 0) continue;
564 
565  // place phi aluminum walls
566  double dphi = aeroGeo.getRingDPhi(iRing);
567 
568  wallName.str("");
569  wallName << "supportWallPhi_" << iRing + 1;
570  G4Box* wall = new G4Box(wallName.str(), (wallR[iRing] - wallR[iRing - 1] - wallThick) / 2. - 1., thick / 2., wallHeight / 2.);
571  G4LogicalVolume* wallLV = new G4LogicalVolume(wall, supportMaterial, string("ARICH.") + wallName.str());
572  double r = (wallR[iRing - 1] + wallThick + wallR[iRing]) / 2.;
573  double zLayer = 0;
574 
575  // loop over layers
576  int iSlot = 1;
577  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
578  double iphi = 0;
579  double layerThick = aeroGeo.getLayerThickness(iLayer);
580 
581  std::stringstream tileName;
582  tileName << "aerogelTile_" << iRing << "_" << iLayer;
583 
584  G4Tubs* tileShape = new G4Tubs(tileName.str(), wallR[iRing - 1] + wallThick + tileGap, wallR[iRing] - tileGap, layerThick / 2.,
585  (tileGap + wallThick / 2.) / wallR[iRing], dphi - (2.*tileGap + wallThick) / wallR[iRing]);
586 
587  G4Material* aeroMaterial = Materials::get(aeroGeo.getLayerMaterial(iLayer));
588  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, aeroMaterial, string("ARICH.") + tileName.str());
589 
590  while (iphi < 2 * M_PI - 0.0001) {
591  G4ThreeVector trans(r * cos(iphi), r * sin(iphi), (thick - imgTubeLen) / 2.);
592  G4RotationMatrix Ra;
593  Ra.rotateZ(iphi);
594 
595  if (iLayer == 1) new G4PVPlacement(G4Transform3D(Ra, trans), wallLV, string("ARICH.") + wallName.str(), aerogelPlaneLV, false,
596  iSlot);
597 
598  G4ThreeVector transTile(0, 0, (thick + layerThick - wallHeight - imgTubeLen) / 2. + zLayer);
599  new G4PVPlacement(G4Transform3D(Ra, transTile), tileLV, string("ARICH.") + tileName.str(), aerogelPlaneLV, false, iSlot);
600  iphi += dphi;
601  iSlot++;
602  }
603  zLayer += layerThick;
604  }
605  }
606 
607  new G4PVPlacement(G4Translate3D(0., 0., -(wallHeight + imgTubeLen) / 2.), supportTubeLV, "ARICH.AerogelSupportPlate",
608  aerogelPlaneLV,
609  false, 1);
610 
611  new G4PVPlacement(G4Translate3D(0., 0., (wallHeight + thick) / 2.), imgTubeLV, "ARICH.AerogelImgPlate", aerogelPlaneLV, false, 1);
612 
613  return aerogelPlaneLV;
614 
615  }
616 
618  {
619 
620  //cout << "GeoARICHCreator::buildAerogelPlaneWithIndividualTilesProp(const ARICHGeometryConfig& detectorGeo)" << endl;
621 
622  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
623 
624  // support plane
625  double rin = aeroGeo.getSupportInnerR();
626  double rout = aeroGeo.getSupportOuterR();
627  double thick = aeroGeo.getSupportThickness();
628  double wallThick = aeroGeo.getWallThickness();
629  // Maximum total (up and down) thickness of the aerogel tiles
630  double maxTotalAerogelThick = aeroGeo.getMaximumTotalTileThickness();
631  //cout<<"maxTotalAerogelThick "<<maxTotalAerogelThick<<endl
632  // <<"wallHeight "<<wallHeight<<endl;
633  // In case of individual thickness of the tiles we need to define compensation
634  // volume with ARICH air. This volume situated between aerogel tile and image plane (imgTube).
635  // Minimum thickness of the compensation volume with ARICH air
636  double compensationARICHairVolumeThick_min = aeroGeo.getCompensationARICHairVolumeThick_min(); // mm
637  // Please note redefinition of wallHeight value
638  double wallHeight = maxTotalAerogelThick + compensationARICHairVolumeThick_min;
639 
640  string supportMat = aeroGeo.getSupportMaterial();
641  G4Material* supportMaterial = Materials::get(supportMat);
642 
643  G4Material* gapMaterial =
644  Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
645  G4Material* imgMaterial = Materials::get("ARICH_Air"); // Air with defined optical properties to propagate cherenkov photons
646 
647  // master volume
648  double imgTubeLen = aeroGeo.getImgTubeThickness(); // if changed, change position of aerogel plane also in main function!!
649  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight + imgTubeLen) / 2., 0, 2 * M_PI);
650  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
651 
652  // support plate
653  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
654  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
655  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
656 
657  // imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
658  G4Tubs* imgTube = new G4Tubs("imgTube", rin, rout, imgTubeLen / 2., 0, 2 * M_PI);
659  G4LogicalVolume* imgTubeLV = new G4LogicalVolume(imgTube, imgMaterial, "ARICH.AerogelImgPlate");
660  imgTubeLV->SetSensitiveDetector(m_sensitiveAero);
661 
662  // read radiuses of aerogel slot aluminum walls
663  std::vector<double> wallR;
664  unsigned nRing = aeroGeo.getNRings();
665  for (unsigned iRing = 1; iRing < nRing + 1; iRing++) {
666  wallR.push_back(aeroGeo.getRingRadius(iRing));
667  }
668 
669  unsigned nLayer = aeroGeo.getNLayers();
670  double tileGap = aeroGeo.getTileGap();
671  G4Transform3D transform = G4Translate3D(0., 0., (thick - imgTubeLen) / 2.);
672 
673  for (unsigned iRing = 0; iRing < nRing; iRing++) {
674 
675  // Aluminum walls between tile rings (r wall)
676  std::stringstream wallName;
677  wallName << "supportWallR_" << iRing + 1;
678  //cout<<"wallName = "<<wallName.str().c_str()<<endl;
679  G4Tubs* supportWall = new G4Tubs(wallName.str().c_str(), wallR[iRing], wallR[iRing] + wallThick, wallHeight / 2., 0, 2 * M_PI);
680  G4LogicalVolume* supportWallLV = new G4LogicalVolume(supportWall, supportMaterial, string("ARICH.") + wallName.str().c_str());
681  new G4PVPlacement(transform, supportWallLV, string("ARICH.") + wallName.str().c_str(), aerogelPlaneLV, false, 0);
683 
684  // There are only 4 rings of aerogel - the first one is only for mechanical support
685  if (iRing == 0) continue;
686 
687  // dphi - distance between centers of two consecutive aluminum walls (diaphragm) in one ring
688  double dphi = aeroGeo.getRingDPhi(iRing);
689 
690  // Aluminum walls (diaphragm) between two neighboring tile in one ring (phi wall)
691  wallName.str("");
692  wallName << "supportWallPhi_" << iRing + 1;
693  G4Box* wall = new G4Box(wallName.str(), (wallR[iRing] - wallR[iRing - 1] - wallThick) / 2. - 1., thick / 2., wallHeight / 2.);
694  G4LogicalVolume* wallLV = new G4LogicalVolume(wall, supportMaterial, string("ARICH.") + wallName.str());
695  double r = (wallR[iRing - 1] + wallThick + wallR[iRing]) / 2.;
697 
698  // loop over layers
699  int icopyNumber = 0;
700  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
701  //int iSlot = 1;
702  double iphi = 0;
703 
704  int iicolumn = 0;
705  // loop over phi (over tile slots from same ring)
706  while (iphi < 2 * M_PI - 0.0001) {
707 
708  // Define layer thickness as -1 in case it will not be defined
709  // below in the code with a appropriate value - Geant4 will trigger error
710  double layerThick = -1.0;
711  double tileUpThick = -1.0;
712  double tileDownThick = -1.0;
713  //cout<<" double layerThick = aeroGeo.getLayerThickness(iLayer) = "<<layerThick<<endl;
714 
715  // Define material and thickness
716  G4Material* aeroMaterial = NULL;
717  int ati_ring = iRing;
718  int ati_column = iicolumn + 1;
719  int ati_layerN = iLayer - 1;
720 
721  //cout<<setw(5)<<ati_layerN<<setw(5)<<ati_ring<<setw(5)<<ati_column<<endl;
722  if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 1) {
723  aeroMaterial = Materials::get(aeroGeo.getTileMaterialName(ati_ring, ati_column, ati_layerN).c_str());
724  } else {
725  B2ERROR("GeoARICHCreator::buildAerogelPlaneWithIndividualTilesProp --> getFullAerogelMaterialDescriptionKey() is wrong");
726  }
727  //cout<<"-----------------"<<endl
728  // <<"iLayer = "<<iLayer<<endl
729  // <<aeroMaterial->GetName()<<endl;
730  //cout<<setw(5)<<ati_layerN<<setw(5)<<ati_ring<<setw(5)<<ati_column<<endl;
731  //aeroMaterial->GetMaterialPropertiesTable()->DumpTable();
732  //zcout<<"ooooooooooooooooo"<<endl;
733  layerThick = aeroGeo.getTileThickness(ati_ring, ati_column, ati_layerN);
734  tileUpThick = aeroGeo.getTileThickness(ati_ring, ati_column, 0);
735  tileDownThick = aeroGeo.getTileThickness(ati_ring, ati_column, 1);
736 
737  // Placement of the aluminum walls (diaphragm) between two neighboring tile in one ring (phi wall)
738  // please note that we place single wall for two aerogel layers
739  G4ThreeVector trans(r * cos(iphi), r * sin(iphi), (thick - imgTubeLen) / 2.);
740  G4RotationMatrix Ra;
741  Ra.rotateZ(iphi);
742  if (iLayer == 1) {
743  new G4PVPlacement(G4Transform3D(Ra, trans), //transformation
744  wallLV, //its logical
745  string("ARICH.") + wallName.str(), //name
746  aerogelPlaneLV, //mother logical
747  false, //always false
748  icopyNumber); //should be set to 0 for the first volume of a given type.
750 
751  // In case of individual thickness
752  // (if aeroGeo.getFullAerogelMaterialDescriptionKey() == 1) or
753  // (if aeroGeo.getFullAerogelMaterialDescriptionKey() == 2)
754  // of the tiles we need to define compensation volume with ARICH air.
755  // This volume situated between aerogel tile and image plane (imgTube).
756  // This volume have same shape as earogel tile but different thickness.
757  // Build Compensation tiles only one time
758  // Compensation tile shape
759  double compTileUpThick = wallHeight - tileUpThick - tileDownThick;
760  std::stringstream compTileName;
761  //compTileName << "aerogelCompTile_" << ati_layerN << "_" << ati_ring << "_" << ati_column;
762  // In the end of the name we have Layer(L), Ring(R), Slot/Column(S) id's
763  // L : 1-2
764  // R : 1-4
765  // S : 1-22 @ R = 1
766  // S : 1-28 @ R = 2
767  // S : 1-34 @ R = 3
768  // S : 1-40 @ R = 4
769  compTileName << "aerogelCompTile_" << iLayer << "_" << ati_ring << "_" << ati_column;
770  //cout<<compTileName.str()<<endl;
771  G4Tubs* compTileShape = new G4Tubs(compTileName.str(), //name
772  wallR[iRing - 1] + wallThick + tileGap, //Rmin
773  wallR[iRing] - tileGap, //Rmax
774  compTileUpThick / 2.0, //Thikness
775  (tileGap + wallThick / 2.0) / wallR[iRing], //phi start
776  dphi - (2.0 * tileGap + wallThick) / wallR[iRing]); //delta phi
777 
778  // Logical volume of the compensation tiles
779  G4LogicalVolume* compTileLV = new G4LogicalVolume(compTileShape, //Its solid
780  imgMaterial, //G4 material
781  string("ARICH.") + compTileName.str()); //name
782 
783  // Placement of the compensation tiles
784  G4ThreeVector transCompTile(0, 0, (thick + wallHeight - compTileUpThick - imgTubeLen) / 2.0);
785  G4RotationMatrix compRa;
786  compRa.rotateZ(iphi);
787  new G4PVPlacement(G4Transform3D(compRa, transCompTile), //transformation
788  compTileLV, //its logical
789  string("ARICH.") + compTileName.str(), //name
790  aerogelPlaneLV, //mother logical
791  false, //always false
792  0); //should be set to 0 for the first volume of a given type.
793  }
794 
795  // Tile shape
796  std::stringstream tileName;
797  tileName << "aerogelTile_" << iLayer << "_" << ati_ring << "_" << ati_column;
798  //cout<<tileName.str()<<endl;
799  G4Tubs* tileShape = new G4Tubs(tileName.str(), //name
800  wallR[iRing - 1] + wallThick + tileGap, //Rmin
801  wallR[iRing] - tileGap, //Rmax
802  layerThick / 2.0, //Thikness
803  (tileGap + wallThick / 2.0) / wallR[iRing], //phi start
804  dphi - (2.0 * tileGap + wallThick) / wallR[iRing]); //delta phi
805 
806  // Logical volume of the aerogel tiles
807  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, //Its solid
808  aeroMaterial, //G4 material
809  string("ARICH.") + tileName.str()); //name
810 
811  // Placement of the aerogel tiles
812  double zLayer = 0.0;
813  if (iLayer == 2)
814  zLayer = tileUpThick;
815  G4ThreeVector transTile(0, 0, (thick + layerThick - wallHeight - imgTubeLen) / 2.0 + zLayer);
816  new G4PVPlacement(G4Transform3D(Ra, transTile), //transformation
817  tileLV, //its logical
818  string("ARICH.") + tileName.str(), //name
819  aerogelPlaneLV, //mother logical
820  false, //always false
821  0); //should be set to 0 for the first volume of a given type.
823 
824  iphi += dphi;
825  //iSlot++;
826  icopyNumber++;
827  iicolumn++;
828  }
829  }
830  }
831 
832  // Placement of the support tube
833  new G4PVPlacement(G4Translate3D(0., 0., -(wallHeight + imgTubeLen) / 2.), //transformation
834  supportTubeLV, //its logical
835  "ARICH.AerogelSupportPlate", //name
836  aerogelPlaneLV, //mother logical
837  false, //always false
838  0); //should be set to 0 for the first volume of a given type.
839 
840  // Placement of the imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
841  new G4PVPlacement(G4Translate3D(0., 0., (wallHeight + thick) / 2.), //transformation
842  imgTubeLV, //its logical
843  "ARICH.AerogelImgPlate", //name
844  aerogelPlaneLV, //mother logical
845  false, //always false
846  0); //should be set to 0 for the first volume of a given type.
847 
848  return aerogelPlaneLV;
849 
850  }
851 
852  G4LogicalVolume* GeoARICHCreator::buildHAPD(const ARICHGeoHAPD& hapdGeo)
853  {
854 
855  // get module materials
856  string wallMat = hapdGeo.getWallMaterial();
857  string winMat = hapdGeo.getWinMaterial();
858  string apdMat = hapdGeo.getAPDMaterial();
859  string fillMat = hapdGeo.getFillMaterial();
860  string febMat = hapdGeo.getFEBMaterial();
861  G4Material* wallMaterial = Materials::get(wallMat);
862  G4Material* windowMaterial = Materials::get(winMat);
863  G4Material* apdMaterial = Materials::get(apdMat);
864  G4Material* fillMaterial = Materials::get(fillMat);
865  G4Material* febMaterial = Materials::get(febMat);
866  G4Material* moduleFill = Materials::get("ARICH_Air");
867 
868  // check that module window material has specified refractive index
869  double wref = getAvgRINDEX(windowMaterial);
870  if (!wref) B2WARNING("Material '" << winMat <<
871  "', required for ARICH photon detector window as no specified refractive index. Continuing, but no photons in ARICH will be detected.");
872 
873  // get module dimensions
874  const double hapdSizeX = hapdGeo.getSizeX();
875  const double hapdSizeY = hapdGeo.getSizeY();
876  const double hapdSizeZ = hapdGeo.getSizeZ();
877  const double wallThick = hapdGeo.getWallThickness();
878  const double winThick = hapdGeo.getWinThickness();
879  const double apdSizeX = hapdGeo.getAPDSizeX();
880  const double apdSizeY = hapdGeo.getAPDSizeY();
881  const double apdSizeZ = hapdGeo.getAPDSizeZ();
882  const double botThick = wallThick;
883  const double modHeight = hapdGeo.getModuleSizeZ();
884 
885  // module master volume
886  G4Box* moduleBox = new G4Box("moduleBox", hapdSizeX / 2., hapdSizeY / 2., modHeight / 2.);
887  G4LogicalVolume* lmoduleBox = new G4LogicalVolume(moduleBox, moduleFill, "ARICH.HAPDModule");
888 
889  // build HAPD box
890  G4Box* hapdBox = new G4Box("hapdBox", hapdSizeX / 2., hapdSizeY / 2., hapdSizeZ / 2.);
891  G4LogicalVolume* lhapdBox = new G4LogicalVolume(hapdBox, fillMaterial, "ARICH.HAPD");
892 
893  // build HAPD walls
894  G4Box* tempBox2 = new G4Box("tempBox2", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick,
895  hapdSizeZ / 2. + 0.1); // Dont't care about "+0.1", needs to be there.
896  G4SubtractionSolid* moduleWall = new G4SubtractionSolid("Box-tempBox", hapdBox, tempBox2);
897  G4LogicalVolume* lmoduleWall = new G4LogicalVolume(moduleWall, wallMaterial, "ARICH.HAPDWall");
898  setColor(*lmoduleWall, "rgb(1.0,0.0,0.0,1.0)");
899  new G4PVPlacement(G4Transform3D(), lmoduleWall, "ARICH.HAPDWall", lhapdBox, false, 1);
900 
901  // build HAPD window
902  G4Box* winBox = new G4Box("winBox", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick, winThick / 2.);
903  G4LogicalVolume* lmoduleWin = new G4LogicalVolume(winBox, windowMaterial, "ARICH.HAPDWindow");
904  setColor(*lmoduleWin, "rgb(0.7,0.7,0.7,1.0)");
905  lmoduleWin->SetSensitiveDetector(m_sensitive);
906  G4Transform3D transform = G4Translate3D(0., 0., (-hapdSizeZ + winThick) / 2.);
907  new G4PVPlacement(transform, lmoduleWin, "ARICH.HAPDWindow", lhapdBox, false, 1);
908 
909  // build module bottom
910  G4Box* botBox = new G4Box("botBox", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick, botThick / 2.);
911  G4LogicalVolume* lmoduleBot = new G4LogicalVolume(botBox, wallMaterial, "ARICH.HAPDBottom");
912  setColor(*lmoduleBot, "rgb(0.0,1.0,0.0,1.0)");
913  G4Transform3D transform1 = G4Translate3D(0., 0., (hapdSizeZ - botThick) / 2.);
914  new G4PVPlacement(transform1, lmoduleBot, "ARICH.HAPDBottom", lhapdBox, false, 1);
915 
916  // build apd
917  G4Box* apdBox = new G4Box("apdBox", apdSizeX / 2., apdSizeY / 2., apdSizeZ / 2.);
918  G4LogicalVolume* lApd = new G4LogicalVolume(apdBox, apdMaterial, "ARICH.HAPDApd");
919  if (m_isBeamBkgStudy) lApd->SetSensitiveDetector(new BkgSensitiveDetector("ARICH", 1));
920 
921  // add APD surface optical properties
922  Materials& materials = Materials::getInstance();
923 
924  G4OpticalSurface* optSurf = materials.createOpticalSurface(hapdGeo.getAPDSurface());
925 
926  new G4LogicalSkinSurface("apdSurface", lApd, optSurf);
927  G4Transform3D transform2 = G4Translate3D(0., 0., (hapdSizeZ - apdSizeZ) / 2. - botThick);
928  new G4PVPlacement(transform2, lApd, "ARICH.HAPDApd", lhapdBox, false, 1);
929 
930  // build FEB
931  double febSizeX = hapdGeo.getFEBSizeX();
932  double febSizeY = hapdGeo.getFEBSizeY();
933  double febSizeZ = hapdGeo.getFEBSizeZ();
934  G4Box* febBox = new G4Box("febBox", febSizeX / 2., febSizeY / 2., febSizeZ / 2.);
935  G4LogicalVolume* lfeb = new G4LogicalVolume(febBox, febMaterial, "ARICH.HAPDFeb");
936  if (m_isBeamBkgStudy) lfeb->SetSensitiveDetector(new BkgSensitiveDetector("ARICH"));
937  setColor(*lfeb, "rgb(0.0,0.6,0.0,1.0)");
938  G4Transform3D transform3 = G4Translate3D(0., 0., (modHeight - febSizeZ) / 2.);
939  new G4PVPlacement(transform3, lfeb, "ARICH.HAPDFeb", lmoduleBox, false, 1);
940  G4Transform3D transform4 = G4Translate3D(0., 0., - (modHeight - hapdSizeZ) / 2.);
941  new G4PVPlacement(transform4, lhapdBox, "ARICH.HAPD", lmoduleBox, false, 1);
942 
943  return lmoduleBox;
944 
945  }
946 
947  G4LogicalVolume* GeoARICHCreator::buildDetectorPlane(const ARICHGeometryConfig& detectorGeo)
948  {
949 
950  const ARICHGeoHAPD& hapdGeo = detectorGeo.getHAPDGeometry();
951  G4LogicalVolume* hapdLV = buildHAPD(hapdGeo);
952 
953  const ARICHGeoDetectorPlane& detGeo = detectorGeo.getDetectorPlane();
954 
955  G4Tubs* detTube = new G4Tubs("detTube", detGeo.getRingR(1) - hapdGeo.getSizeX() * 1.4 / 2.,
956  detGeo.getRingR(detGeo.getNRings()) + hapdGeo.getSizeX() * 1.4 / 2., hapdGeo.getModuleSizeZ() / 2., 0, 2 * M_PI);
957  G4LogicalVolume* detPlaneLV = new G4LogicalVolume(detTube, Materials::get("ARICH_Air"), "ARICH.detectorPlane");
958 
959  unsigned nSlots = detGeo.getNSlots();
960 
961  for (unsigned iSlot = 1; iSlot < nSlots + 1; iSlot++) {
962  if (!m_modInfo->isInstalled(iSlot)) continue;
963  double r = detGeo.getSlotR(iSlot);
964  double phi = detGeo.getSlotPhi(iSlot);
965  G4ThreeVector trans(r * cos(phi), r * sin(phi), 0);
966  G4RotationMatrix Ra;
967  Ra.rotateZ(phi);
968  G4ThreeVector trans1(r * cos(phi), r * sin(phi), 0.0);
969  new G4PVPlacement(G4Transform3D(Ra, trans1), hapdLV, "ARICH.HAPDModule", detPlaneLV, false, iSlot);
970  }
971 
972  return detPlaneLV;
973 
974  }
975 
976  G4LogicalVolume* GeoARICHCreator::buildMerger(const ARICHGeoMerger& mergerGeo)
977  {
978 
979  // This is a screw hole on the merger board.
980  // This high precision of the geometry description is needed
981  // only for the correct placement of the merger cooling bodies.
982  // Volume to subtract (screw hole on the merger board)
983  double screwholeR = mergerGeo.getMergerPCBscrewholeR() * mm;
984  double screwholedY = mergerGeo.getMergerPCBscrewholePosdY() * mm;
985  double screwholedX1 = mergerGeo.getMergerPCBscrewholePosdX1() * mm;
986  double screwholedX2 = mergerGeo.getMergerPCBscrewholePosdX2() * mm;
987 
988  G4VSolid* screwHoleTubeSubtracted_solid = new G4Tubs("screwHoleTubeSubtracted_solid",
989  0.0,
990  screwholeR,
991  mergerGeo.getThickness() * mm / 2.0,
992  0, 360.0 * deg);
993 
994  // Volume to add (merger box)
995  G4Box* merger_solid = new G4Box("merger_solid",
996  mergerGeo.getSizeW() * mm / 2.0,
997  mergerGeo.getSizeL() * mm / 2.0,
998  mergerGeo.getThickness() * mm / 2.0);
999 
1000  G4RotationMatrix Ra_sub;
1001  G4ThreeVector Ta_sub;
1002  G4Transform3D Tr_sub;
1003  Ta_sub.setX(-mergerGeo.getSizeW() * mm / 2.0 + screwholedX1);
1004  Ta_sub.setY(mergerGeo.getSizeL() * mm / 2.0 - screwholedY);
1005  Ta_sub.setZ(0.0);
1006  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1007  G4SubtractionSolid* substraction_solid = new G4SubtractionSolid("substraction_solid", merger_solid, screwHoleTubeSubtracted_solid,
1008  Tr_sub);
1009  Ta_sub.setX(mergerGeo.getSizeW() * mm / 2.0 - screwholedX2);
1010  Ta_sub.setY(mergerGeo.getSizeL() * mm / 2.0 - screwholedY);
1011  Ta_sub.setZ(0.0);
1012  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1013  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, screwHoleTubeSubtracted_solid, Tr_sub);
1014  Ta_sub.setX(mergerGeo.getSizeW() * mm / 2.0 - screwholedX2);
1015  Ta_sub.setY(-mergerGeo.getSizeL() * mm / 2.0 + screwholedY);
1016  Ta_sub.setZ(0.0);
1017  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1018  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, screwHoleTubeSubtracted_solid, Tr_sub);
1019  Ta_sub.setX(-mergerGeo.getSizeW() * mm / 2.0 + screwholedX1);
1020  Ta_sub.setY(-mergerGeo.getSizeL() * mm / 2.0 + screwholedY);
1021  Ta_sub.setZ(0.0);
1022  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1023  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, screwHoleTubeSubtracted_solid, Tr_sub);
1024 
1025  return new G4LogicalVolume(substraction_solid, Materials::get(mergerGeo.getMergerPCBMaterialName()), "ARICH.mergerPCB");
1026  }
1027 
1028 
1029  G4LogicalVolume* GeoARICHCreator::buildMergerCooling(unsigned iType)
1030  {
1031 
1032  if (!m_mergerCooling) {
1033  B2WARNING("ARICH geometry: no data available for merger " << iType << " cooling body geometry. Cooling body will not be placed.");
1034  return NULL;
1035  }
1036 
1037  std::stringstream shpName;
1038  shpName << "TessellatedSolid_" << + iType;
1039 
1040  G4TessellatedSolid* volume_solid = new G4TessellatedSolid(shpName.str().c_str());
1041 
1042  G4ThreeVector point_1;
1043  G4ThreeVector point_2;
1044  G4ThreeVector point_3;
1045 
1046  tessellatedSolidStr mergerCoolingStr = m_mergerCooling->getMergerCoolingBodiesInfo(iType);
1047 
1048  if (mergerCoolingStr.nCells == 0) {
1049  B2WARNING("ARICH geometry: no data available for merger " << iType << " cooling body geometry. Cooling body will not be placed.");
1050  return NULL;
1051  }
1052 
1053  for (unsigned int i = 0; i < mergerCoolingStr.nCells; i++) {
1054  //
1055  point_1.setX(mergerCoolingStr.posV1[0][i]);
1056  point_1.setY(mergerCoolingStr.posV1[1][i]);
1057  point_1.setZ(mergerCoolingStr.posV1[2][i]);
1058  //
1059  point_2.setX(mergerCoolingStr.posV2[0][i]);
1060  point_2.setY(mergerCoolingStr.posV2[1][i]);
1061  point_2.setZ(mergerCoolingStr.posV2[2][i]);
1062  //
1063  point_3.setX(mergerCoolingStr.posV3[0][i]);
1064  point_3.setY(mergerCoolingStr.posV3[1][i]);
1065  point_3.setZ(mergerCoolingStr.posV3[2][i]);
1066  //
1067  G4TriangularFacet* facet = new G4TriangularFacet(point_1, point_2, point_3, ABSOLUTE);
1068  volume_solid->AddFacet((G4VFacet*) facet);
1069  delete facet;
1070  }
1071 
1072  volume_solid->SetSolidClosed(true);
1073  std::stringstream volName;
1074  volName << "ARICH.mergerCooling_" << + iType;
1075  G4LogicalVolume* volume_logical = new G4LogicalVolume(volume_solid,
1076  Materials::get(m_mergerCooling->getMergerCoolingBodiesMaterialName()), volName.str().c_str());
1077 
1078  setColor(*volume_logical, "rgb(0.6,0.0,0.2,1.0)"); //From the top (farther from IP, downstream)
1079 
1080  return volume_logical;
1081  }
1082 
1083  G4LogicalVolume* GeoARICHCreator::buildMergerEnvelope(const ARICHGeoMerger& mergerGeo, int iType)
1084  {
1085  // Volume to add single merger and merger cooling body envelope box
1086  G4Box* singlemergerenvelope_solid = new G4Box("singlemergerenvelope_solid",
1087  mergerGeo.getSingleMergerEnvelopeSizeW() * mm / 2.0,
1088  mergerGeo.getSingleMergerEnvelopeSizeL() * mm / 2.0,
1089  mergerGeo.getSingleMergerEnvelopeThickness() * mm / 2.0);
1090  std::stringstream volName;
1091  volName << "ARICH.singleMergerEnvelope_" << + iType;
1092  return new G4LogicalVolume(singlemergerenvelope_solid, Materials::get("ARICH_Air"), volName.str().c_str());
1093  }
1094 
1096  {
1097 
1098  const ARICHGeoMerger& mergerGeo = detectorGeo.getMergerGeometry();
1099 
1100  if (mergerGeo.getSingleMergerEnvelopeSizeW() < 1e-9) {
1101  B2WARNING("GeoARICHCreator: Merger and merger cooling geometry will not be build as it is not availible in geometry configuration (ARICHGeometryConfig with ClasDef>4 is needed).");
1102  return NULL;
1103  }
1104 
1105  G4Tubs* envelope_solid = new G4Tubs("envelope_solid", mergerGeo.getEnvelopeInnerRadius() * mm,
1106  mergerGeo.getEnvelopeOuterRadius() * mm, mergerGeo.getEnvelopeThickness() * mm / 2.0, 0.0, 2.0 * M_PI);
1107 
1108  G4LogicalVolume* envelope_logical = new G4LogicalVolume(envelope_solid, Materials::get("ARICH_Air"), "ARICH.mergerEnvelope");
1109 
1110 
1111  G4LogicalVolume* merger_logical = buildMerger(mergerGeo);
1112  G4LogicalVolume* mergerCooling_logical[12] = {NULL};
1113  G4LogicalVolume* mergerEnvelope_logical[12] = {NULL};
1114  G4ThreeVector TaPCB(mergerGeo.getSingleMergeEnvelopePosition().X() * mm, mergerGeo.getSingleMergeEnvelopePosition().Y() * mm,
1115  mergerGeo.getSingleMergeEnvelopePosition().Z() * mm);
1116  G4RotationMatrix RaPCB;
1117  G4ThreeVector TaMergerCooling(mergerGeo.getSingleMergeEnvelopePosition().X() * mm,
1118  mergerGeo.getSingleMergeEnvelopePosition().Y() * mm, -mergerGeo.getSingleMergeEnvelopePosition().Z() * mm);
1119  G4RotationMatrix RaMergerCooling;
1120  RaMergerCooling.rotateY(180 * deg);
1121  RaMergerCooling.rotateZ(-90 * deg);
1122 
1123  // build 12 different merger+cooling body volumes
1124  for (int iType = 1; iType < 13; iType++) {
1125  mergerCooling_logical[iType - 1] = buildMergerCooling(iType);
1126  mergerEnvelope_logical[iType - 1 ] = buildMergerEnvelope(mergerGeo, iType); //Single merger and merger cooling body envelope box
1127  setColor(*mergerEnvelope_logical[iType - 1], "rgb(0.0,0.0,1.0,1.0)");
1128 
1129  new G4PVPlacement(G4Transform3D(RaPCB, TaPCB), //Transformation
1130  merger_logical, //its logical volume
1131  "ARICH.mergerPCB", //its name
1132  mergerEnvelope_logical[iType - 1], //its mother volume
1133  false, //no boolean operation
1134  iType); //copy number
1135 
1136  if (mergerCooling_logical[iType - 1] == NULL) continue;
1137 
1138  new G4PVPlacement(G4Transform3D(RaMergerCooling, TaMergerCooling), //Transformation
1139  mergerCooling_logical[iType - 1], //its logical volume
1140  "ARICH.mergerCooling", //its name
1141  mergerEnvelope_logical[iType - 1], //its mother volume
1142  false, //no boolean operation
1143  iType); //copy number
1144  }
1145 
1146  // place all 72 merger+cooling body packages
1147  for (unsigned iSlot = 0; iSlot < mergerGeo.getMergerSlotID().size(); iSlot++) {
1148 
1149  int type = 1; // if no merger cooling is available...
1150  if (m_mergerCooling) type = (int)m_mergerCooling->getMergerCoolingPositionID().at(iSlot);
1151 
1152  //Placement of the single volume envelope
1153  G4ThreeVector Ta(mergerGeo.getMergerPosR().at(iSlot) * mm * cos(mergerGeo.getMergerAngle().at(iSlot) * deg),
1154  mergerGeo.getMergerPosR().at(iSlot) * mm * sin(mergerGeo.getMergerAngle().at(iSlot) * deg),
1155  mergerGeo.getSingleMergerenvelopeDeltaZ().at(iSlot) * mm);
1156  G4RotationMatrix Ra;
1157  Ra.rotateZ(mergerGeo.getMergerAngle().at(iSlot) * deg + mergerGeo.getMergerOrientation().at(iSlot) * deg);
1158  new G4PVPlacement(G4Transform3D(Ra, Ta), //Transformation
1159  mergerEnvelope_logical[type - 1], //its logical volume
1160  "ARICH.singleMergerEnvelope", //its name
1161  envelope_logical, //its mother volume
1162  false, //no boolean operation
1163  iSlot); //copy number
1164  }
1165 
1166  return envelope_logical;
1167 
1168  }
1169 
1170  G4LogicalVolume* GeoARICHCreator::buildCables(const ARICHGeoCablesEnvelope& cablesGeo)
1171  {
1172 
1173  G4Tubs* cablesEnvelope_solid = new G4Tubs("cablesEnvelope_solid",
1174  cablesGeo.getEnvelopeInnerRadius(),
1175  cablesGeo.getEnvelopeOuterRadius(),
1176  cablesGeo.getEnvelopeThickness() / 2.0,
1177  0.0,
1178  2.0 * M_PI);
1179  G4LogicalVolume* cablesEnvelope_logical = new G4LogicalVolume(cablesEnvelope_solid,
1181  "ARICH.cablesEnvelope");
1182 
1183  return cablesEnvelope_logical;
1184 
1185  }
1186 
1187  G4LogicalVolume* GeoARICHCreator::buildFEBCoolingBody(const ARICHGeoFEBCooling& coolingv2Geo)
1188  {
1189 
1190  //FEB aluminum cooling envelope for single object
1191  double feb_alcooling_singleObjectEnvelope_sizeX = (2 * coolingv2Geo.getSmallSquareSize() + coolingv2Geo.getBigSquareSize() + 2.0) *
1192  mm;
1193  double feb_alcooling_singleObjectEnvelope_sizeY = feb_alcooling_singleObjectEnvelope_sizeX * mm;
1194  double feb_alcooling_singleObjectEnvelope_sizeZ = (coolingv2Geo.getSmallSquareThickness() + coolingv2Geo.getRectangleThickness()) *
1195  mm;
1196 
1197  double feb_alcooling_box1_sizeX = coolingv2Geo.getSmallSquareSize() * mm;
1198  double feb_alcooling_box1_sizeY = feb_alcooling_box1_sizeX;
1199  double feb_alcooling_box1_sizeZ = coolingv2Geo.getSmallSquareThickness() * mm;
1200 
1201  double feb_alcooling_box2_sizeX = coolingv2Geo.getBigSquareSize() * mm;
1202  double feb_alcooling_box2_sizeY = feb_alcooling_box2_sizeX;
1203  double feb_alcooling_box2_sizeZ = coolingv2Geo.getBigSquareThickness() * mm;
1204 
1205  double feb_alcooling_box3_sizeX = coolingv2Geo.getRectangleW() * mm;
1206  double feb_alcooling_box3_sizeY = coolingv2Geo.getRectangleL() * mm;
1207  double feb_alcooling_box3_sizeZ = coolingv2Geo.getRectangleThickness() * mm;
1208 
1209  double feb_alcooling_box1_X0 = feb_alcooling_box2_sizeX / 2.0 + feb_alcooling_box1_sizeX / 2.0;
1210  double feb_alcooling_box1_Y0 = feb_alcooling_box2_sizeY / 2.0 + feb_alcooling_box1_sizeY / 2.0;
1211  double feb_alcooling_box1_Z0 = 0.0 * mm;
1212 
1213  //double feb_alcooling_box2_X0 = 0.0*mm;
1214  //double feb_alcooling_box2_Y0 = 0.0*mm;
1215  //double feb_alcooling_box2_Z0 = 0.0*mm;
1216 
1217  double feb_alcooling_box3_X0 = coolingv2Geo.getRectangleDistanceFromCenter() / sqrt(2.0) * mm;
1218  double feb_alcooling_box3_Y0 = feb_alcooling_box3_X0;
1219  double feb_alcooling_box3_Z0 = feb_alcooling_box1_sizeZ / 2.0 + feb_alcooling_box3_sizeZ / 2.0;
1220  double feb_alcooling_box3_angle = 45.0 * deg;
1221 
1222  G4RotationMatrix Ra;
1223  G4ThreeVector Ta;
1224  G4Transform3D Tr;
1225 
1226  //
1227  // Define single FEB aluminum cooling envelope
1228  //
1229  G4VSolid* feb_alcoolingEnvelope_solid = new G4Box("feb_alcoolingEnvelope_solid",
1230  feb_alcooling_singleObjectEnvelope_sizeX / 2.0,
1231  feb_alcooling_singleObjectEnvelope_sizeY / 2.0,
1232  feb_alcooling_singleObjectEnvelope_sizeZ / 2.0);
1233  G4LogicalVolume* feb_alcoolingEnvelope_logical = new G4LogicalVolume(feb_alcoolingEnvelope_solid, Materials::get("Air"),
1234  "feb_alcoolingEnvelope_logical");
1235 
1236  G4VSolid* feb_alcooling_box1_solid = new G4Box("feb_alcooling_box1_solid", feb_alcooling_box1_sizeX / 2.0,
1237  feb_alcooling_box1_sizeY / 2.0, feb_alcooling_box1_sizeZ / 2.0);
1238  G4VSolid* feb_alcooling_box2_solid = new G4Box("feb_alcooling_box2_solid", feb_alcooling_box2_sizeX / 2.0,
1239  feb_alcooling_box2_sizeY / 2.0, feb_alcooling_box2_sizeZ / 2.0);
1240  G4VSolid* feb_alcooling_box3_solid = new G4Box("feb_alcooling_box3_solid", feb_alcooling_box3_sizeX / 2.0,
1241  feb_alcooling_box3_sizeY / 2.0, feb_alcooling_box3_sizeZ / 2.0);
1242 
1243  //
1244  //Box 1 A
1245  //
1246  Ta.setX(feb_alcooling_box1_X0);
1247  Ta.setY(feb_alcooling_box1_Y0);
1248  Ta.setZ(feb_alcooling_box1_Z0);
1249  Tr = G4Transform3D(Ra, Ta);
1250  G4UnionSolid* feb_alcooling_assembly01_solid = new G4UnionSolid("feb_alcooling_assembly01_solid", feb_alcooling_box2_solid,
1251  feb_alcooling_box1_solid, Tr);
1252  //
1253  //Box 1 B
1254  //
1255  Ta.setX(-feb_alcooling_box1_X0);
1256  Ta.setY(-feb_alcooling_box1_Y0);
1257  Ta.setZ(feb_alcooling_box1_Z0);
1258  Tr = G4Transform3D(Ra, Ta);
1259  G4UnionSolid* feb_alcooling_assembly02_solid = new G4UnionSolid("feb_alcooling_assembly02_solid", feb_alcooling_assembly01_solid,
1260  feb_alcooling_box1_solid, Tr);
1261  //
1262  //Box 3 A
1263  //
1264  Ta.setX(feb_alcooling_box3_X0);
1265  Ta.setY(feb_alcooling_box3_Y0);
1266  Ta.setZ(feb_alcooling_box3_Z0);
1267  Ra.rotateZ(-feb_alcooling_box3_angle);
1268  Tr = G4Transform3D(Ra, Ta);
1269  G4UnionSolid* feb_alcooling_assembly03_solid = new G4UnionSolid("feb_alcooling_assembly03_solid", feb_alcooling_assembly02_solid,
1270  feb_alcooling_box3_solid, Tr);
1271  Ra.rotateZ(feb_alcooling_box3_angle);
1272  //
1273  //Box 3 B
1274  //
1275  Ta.setX(-feb_alcooling_box3_X0);
1276  Ta.setY(-feb_alcooling_box3_Y0);
1277  Ta.setZ(feb_alcooling_box3_Z0);
1278  Ra.rotateZ(-feb_alcooling_box3_angle);
1279  Tr = G4Transform3D(Ra, Ta);
1280  G4UnionSolid* feb_alcooling_assembly_solid = new G4UnionSolid("feb_alcooling_assembly_solid", feb_alcooling_assembly03_solid,
1281  feb_alcooling_box3_solid, Tr);
1282  Ra.rotateZ(feb_alcooling_box3_angle);
1283 
1284  G4LogicalVolume* feb_alcooling_assembly_logical = new G4LogicalVolume(feb_alcooling_assembly_solid, Materials::get("Al"),
1285  "feb_alcooling_assembly_logical");
1286  Ta.setX(0.0);
1287  Ta.setY(0.0);
1288  Ta.setZ(-feb_alcooling_box3_sizeZ / 2.0);
1289 
1290  Tr = G4Transform3D(Ra, Ta);
1291  new G4PVPlacement(Tr, //Transformation
1292  feb_alcooling_assembly_logical, //its logical volume
1293  "feb_alcooling_assembly", //its name
1294  feb_alcoolingEnvelope_logical, //its mother volume
1295  false, //no boolean operation
1296  0); //copy number
1297 
1298  return feb_alcoolingEnvelope_logical;
1299  }
1300 
1301  G4LogicalVolume* GeoARICHCreator::buildCoolingTube(const unsigned i_volumeID, const ARICHGeoCooling& coolingGeo)
1302  {
1303 
1304  B2ASSERT("ARICH cooling geometry ID (G4Tube) is wrong : coolingGeo.getCoolingGeometryID.at(i_volumeID) != 1",
1305  coolingGeo.getCoolingGeometryID().at(i_volumeID) == 1);
1306  G4Tubs* coolingTube_solid = new G4Tubs("coolingTube_solid",
1307  coolingGeo.getRmin() * mm,
1308  coolingGeo.getRmax() * mm,
1309  coolingGeo.getCoolingL().at(i_volumeID) * mm / 2.0,
1310  0.0,
1311  2.0 * M_PI);
1312  return new G4LogicalVolume(coolingTube_solid, Materials::get(coolingGeo.getCoolingPipeMaterialName()), "ARICH.coolingTube");
1313 
1314  }
1315 
1316  G4LogicalVolume* GeoARICHCreator::buildCoolingTorus(const unsigned i_volumeID, const ARICHGeoCooling& coolingGeo)
1317  {
1318 
1319  B2ASSERT("ARICH cooling geometry ID (G4Torus) is wrong : coolingGeo.getCoolingGeometryID.at(i_volumeID) != 2",
1320  coolingGeo.getCoolingGeometryID().at(i_volumeID) == 2);
1321 
1322  double pSPhi = coolingGeo.getCoolingPosPhi().at(i_volumeID) * deg - coolingGeo.getCoolingL().at(
1323  i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID) / 2.0;
1324  //double pDPhi = coolingGeo.getCoolingPosPhi().at(i_volumeID)*deg + coolingGeo.getCoolingL().at(i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID) / 2.0;
1325  double pDPhi = coolingGeo.getCoolingL().at(i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID);
1326 
1327  G4Torus* coolingTorus_solid = new G4Torus("coolingTorus_solid", // name
1328  coolingGeo.getRmin(), // pRmin
1329  coolingGeo.getRmax(), // pRmax
1330  coolingGeo.getCoolingPosR().at(i_volumeID), // pRtor
1331  pSPhi, // pSPhi
1332  pDPhi); // pDPhi
1333 
1334  return new G4LogicalVolume(coolingTorus_solid, Materials::get(coolingGeo.getCoolingPipeMaterialName()), "ARICH.coolingTorus");
1335 
1336  }
1337 
1339  {
1340 
1341  G4Tubs* coolingEnvelope_solid = new G4Tubs("coolingEnvelope_solid",
1342  coolingGeo.getEnvelopeInnerRadius(),
1343  coolingGeo.getEnvelopeOuterRadius(),
1344  coolingGeo.getEnvelopeThickness() / 2.0,
1345  0.0,
1346  2.0 * M_PI);
1347  G4LogicalVolume* coolingEnvelope_logical = new G4LogicalVolume(coolingEnvelope_solid,
1348  Materials::get("Air"),
1349  "ARICH.coolingEnvelope");
1350 
1351  unsigned nComponents = coolingGeo.getCoolingGeometryID().size();
1352 
1353  for (unsigned i = 0; i < nComponents; i++) {
1354  double r = coolingGeo.getCoolingPosR().at(i) * mm;
1355  double phi = coolingGeo.getCoolingPosPhi().at(i) * deg;
1356  G4ThreeVector Ta;
1357  G4RotationMatrix Ra;
1358  G4LogicalVolume* coolingComponentLV;
1359  if (coolingGeo.getCoolingGeometryID().at(i) == 1) {
1360  //<!-- 1 -> G4Tubs --> build a tube
1361  Ta.set(r * cos(phi), r * sin(phi), 0);
1362  //First need to rotate around y axis
1363  //to make the tube parallel with x axis
1364  Ra.rotateY(90.0 * deg);
1365  Ra.rotateZ(coolingGeo.getCoolinRotationAngle().at(i) * deg);
1366  coolingComponentLV = buildCoolingTube(i, coolingGeo);
1367  } else if (coolingGeo.getCoolingGeometryID().at(i) == 2) {
1368  //<!-- 2 -> G4Torus --> build a torus
1369  coolingComponentLV = buildCoolingTorus(i, coolingGeo);
1370  } else {
1371  B2FATAL("ARICH cooling geometry component ID is wrong");
1372  }
1373  new G4PVPlacement(G4Transform3D(Ra, Ta), //Transformation
1374  coolingComponentLV, //its logical volume
1375  "ARICH.cooling", //its name
1376  coolingEnvelope_logical, //its mother volume
1377  false, //no boolean operation
1378  i); //copy number
1379  }
1380 
1381  return coolingEnvelope_logical;
1382 
1383  }
1384 
1385  G4LogicalVolume* GeoARICHCreator::buildCoolingTestPlate(const ARICHGeoCooling& coolingGeo)
1386  {
1387 
1388  G4Box* coolingTestPlateEnvelop_solid = new G4Box("coolingTestPlateEnvelop_solid",
1389  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1390  coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 * mm,
1391  coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 * mm);
1392  G4LogicalVolume* coolingTestPlateEnvelop_logical = new G4LogicalVolume(coolingTestPlateEnvelop_solid, Materials::get("Air"),
1393  "ARICH.coolingTestPlateEnvelop");
1394 
1395  G4Box* coolingTestPlate_solid = new G4Box("coolingTestPlate_solid",
1396  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1397  coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 * mm,
1398  coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 * mm);
1399 
1400  // Volume to subtract
1401  G4VSolid* coldTubeSubtracted_solid = new G4Tubs("coldTubeSubtracted_solid",
1402  0.0,
1403  coolingGeo.getColdTubeSubtractedR() * mm,
1404  (coolingGeo.getCoolingTestPlateslengths().X() + 1) / 2.0 * mm,
1405  0, 360.0 * deg);
1406 
1407  // Volume to add (cold tube)
1408  G4VSolid* coldTube_solid = new G4Tubs("coldTube_solid",
1409  (coolingGeo.getColdTubeR() - coolingGeo.getColdTubeWallThickness()) * mm,
1410  coolingGeo.getColdTubeR() * mm,
1411  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1412  0, 360.0 * deg);
1413  G4LogicalVolume* coldTube_logical = new G4LogicalVolume(coldTube_solid, Materials::get(coolingGeo.getColdTubeMaterialName()),
1414  "ARICH.coldTube");
1415 
1416  G4RotationMatrix Ra_sub;
1417  G4ThreeVector Ta_sub;
1418  G4Transform3D Tr_sub;
1419  Ta_sub.setX(0.0);
1420  Ta_sub.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing()) * mm);
1421  Ta_sub.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1422  Ra_sub.rotateY(90.0 * deg);
1423  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1424  G4SubtractionSolid* substraction_solid = new G4SubtractionSolid("substraction_solid", coolingTestPlate_solid,
1425  coldTubeSubtracted_solid, Tr_sub);
1426  for (int i = 1; i < coolingGeo.getColdTubeNumber(); i++) {
1427  Ta_sub.setX(0.0);
1428  Ta_sub.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing() - coolingGeo.getColdTubeSpacing()
1429  * 2 * i) * mm);
1430  Ta_sub.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1431  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1432  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, coldTubeSubtracted_solid, Tr_sub);
1433  }
1434 
1435  G4LogicalVolume* coolingTestPlate_logical = new G4LogicalVolume(substraction_solid,
1436  Materials::get(coolingGeo.getCoolingTestPlateMaterialName()), "ARICH.coolingTestPlate");
1437 
1438  new G4PVPlacement(G4Transform3D(), //Transformation
1439  coolingTestPlate_logical, //its logical volume
1440  "ARICH.coolingTestPlate", //its name
1441  coolingTestPlateEnvelop_logical, //its mother volume
1442  false, //no boolean operation
1443  0); //copy number
1444  //Add cold tubes
1445  G4RotationMatrix Ra;
1446  G4ThreeVector Ta;
1447  G4Transform3D Tr;
1448  Ra.rotateY(90.0 * deg);
1449  for (int i = 0; i < coolingGeo.getColdTubeNumber(); i++) {
1450  Ta.setX(0.0);
1451  Ta.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing() - coolingGeo.getColdTubeSpacing() *
1452  2 * i) * mm);
1453  Ta.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1454  Tr = G4Transform3D(Ra, Ta);
1455  new G4PVPlacement(Tr, //Transformation
1456  coldTube_logical, //its logical volume
1457  "ARICH.coldTube", //its name
1458  coolingTestPlateEnvelop_logical, //its mother volume
1459  false, //no boolean operation
1460  0); //copy number
1461  }
1462 
1463  return coolingTestPlateEnvelop_logical;
1464 
1465  }
1466 
1468  {
1469 
1470  const ARICHGeoDetectorPlane& detGeo = detectorGeo.getDetectorPlane();
1471 
1472  G4Tubs* supportTube = new G4Tubs("supportTube", detGeo.getSupportInnerR(), detGeo.getSupportOuterR(),
1473  (detGeo.getSupportThickness() + detGeo.getSupportBackWallHeight()) / 2., 0, 2 * M_PI);
1474  G4Material* supportMaterial = Materials::get(detGeo.getSupportMaterial());
1475  G4LogicalVolume* detSupportLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.detectorSupportPlate");
1476 
1477  G4Tubs* supportPlate = new G4Tubs("supportPlate", detGeo.getSupportInnerR(), detGeo.getSupportOuterR(),
1478  detGeo.getSupportThickness() / 2., 0, 2 * M_PI);
1479 
1480  G4Box* hole = new G4Box("hole", detGeo.getModuleHoleSize() / 2., detGeo.getModuleHoleSize() / 2.,
1481  detGeo.getSupportThickness() / 2.); // +1 for thickness for subtraction solid
1482  G4LogicalVolume* holeLV = new G4LogicalVolume(hole, Materials::get("Air"), "ARICH.detectorSupportHole");
1483 
1484  int nRings = detGeo.getNRings();
1485  std::vector<G4LogicalVolume*> hapdBackRadialWallLV;
1486  double backWallThick = detGeo.getSupportBackWallThickness();
1487  double backWallHeight = detGeo.getSupportBackWallHeight();
1488 
1489  G4LogicalVolume* hapdSupportPlateLV = new G4LogicalVolume(supportPlate, supportMaterial, "hapdSupport");
1490 
1491  std::vector<double> wallR;
1492  wallR.assign(nRings + 1, 0);
1493  std::vector<double> thickR;
1494  thickR.assign(nRings + 1, 0);
1495 
1496  for (int i = 1; i < nRings; i++) {
1497  double rm1 = detGeo.getRingR(i);
1498  double rp1 = detGeo.getRingR(i + 1);
1499  wallR[i] = (rp1 + rm1) / 2.;
1500  if (i == 1) {
1501  wallR[0] = rm1 - (rp1 - rm1) / 2.;
1502  }
1503  if (i == nRings - 1) {
1504  wallR[i + 1] = rp1 + (rp1 - rm1) / 2.;
1505  }
1506  }
1507 
1508  for (int i = 0; i < nRings + 1; i++) {
1509  std::stringstream ringName1;
1510  ringName1 << "backWall_" << i;
1511  thickR[i] = backWallThick;
1512  if (i == nRings) {
1513  thickR[i] = 2.*backWallThick;
1514  wallR[i] = detGeo.getSupportOuterR() - thickR[i] / 2.;
1515  }
1516  G4Tubs* backTube = new G4Tubs("hapdBackRing", wallR[i] - thickR[i] / 2., wallR[i] + thickR[i] / 2., backWallHeight / 2., 0,
1517  2 * M_PI);
1518  G4LogicalVolume* hapdBackTubeLV = new G4LogicalVolume(backTube, supportMaterial, "backTube");
1519  G4Transform3D transform3 = G4Translate3D(0., 0., detGeo.getSupportThickness() / 2.);
1520  new G4PVPlacement(transform3, hapdBackTubeLV, "backTube", detSupportLV, false, 1);
1521  if (i == 0) continue;
1522 
1523  G4Box* backRadial = new G4Box("backRadialBox", (wallR[i] - wallR[i - 1] - thickR[i] / 2. - thickR[i - 1] / 2.) / 2. - 1.,
1524  backWallThick / 2., backWallHeight / 2.);
1525  hapdBackRadialWallLV.push_back(new G4LogicalVolume(backRadial, supportMaterial, ringName1.str().c_str()));
1526  }
1527 
1528  G4SubtractionSolid* substraction = NULL;
1529  unsigned nSlots = detGeo.getNSlots();
1530 
1531  if (nSlots > detectorGeo.getFEBCoolingGeometry().getFebcoolingv2GeometryID().size()) {
1532  B2WARNING("GeoARICHCreator: No FEB colling body geometry available so they will not be placed (ARICHGeometryConfig with ClasDef>4 is needed).");
1533  return detSupportLV;
1534  }
1535 
1536  // drill holes in support plate
1537  // add FEB cooling bodies
1538  for (unsigned iSlot = 1; iSlot < nSlots + 1; iSlot++) {
1539  unsigned iRing = detGeo.getSlotRing(iSlot);
1540  double r = (wallR[iRing] + wallR[iRing - 1]) / 2. - (thickR[iRing] - thickR[iRing - 1]) / 2.;
1541 
1542  double phi = detGeo.getSlotPhi(iSlot);
1543  G4ThreeVector trans(r * cos(phi), r * sin(phi), 0);
1544  G4RotationMatrix Ra;
1545  Ra.rotateZ(phi);
1546 
1547  new G4PVPlacement(G4Transform3D(Ra, trans), holeLV, "hole", hapdSupportPlateLV, false, iSlot);
1548  if (substraction) substraction = new G4SubtractionSolid("Box+CylinderMoved", substraction, hole, G4Transform3D(Ra, trans));
1549  else substraction = new G4SubtractionSolid("Box+CylinderMoved", supportPlate, hole, G4Transform3D(Ra, trans));
1550 
1551  phi = phi + detGeo.getRingDPhi(iRing) / 2.;
1552  G4ThreeVector transBack(r * cos(phi), r * sin(phi), detGeo.getSupportThickness() / 2.);
1553  G4RotationMatrix RaBack;
1554  RaBack.rotateZ(phi);
1555  new G4PVPlacement(G4Transform3D(RaBack, transBack), hapdBackRadialWallLV[iRing - 1], "hapdBack", detSupportLV, false, iSlot);
1556 
1557  // add FEB cooling bodies
1558  G4ThreeVector febCoolingTa;
1559  G4Transform3D febCoolingTr;
1560  febCoolingTa.setX(r * cos(detGeo.getSlotPhi(iSlot)));
1561  febCoolingTa.setY(r * sin(detGeo.getSlotPhi(iSlot)));
1562 
1563  double supportTube_envelope_dZ = (detGeo.getSupportThickness() + detGeo.getSupportBackWallHeight());
1564  double febCooling_envelope_dZ = (detectorGeo.getFEBCoolingGeometry().getBigSquareThickness() +
1566  double febCooling_envelope_Z0 = -supportTube_envelope_dZ / 2.0 + febCooling_envelope_dZ / 2.0 + detGeo.getSupportThickness();
1567  febCoolingTa.setZ(febCooling_envelope_Z0);
1568 
1569  int febcoolingv2GeometryID = detectorGeo.getFEBCoolingGeometry().getFebcoolingv2GeometryID().at(iSlot - 1);
1570 
1571  if (febcoolingv2GeometryID == 2) Ra.rotateZ(90.0 * deg);
1572 
1573  febCoolingTr = G4Transform3D(Ra, febCoolingTa);
1574 
1575  if (febcoolingv2GeometryID != 0) {
1576  G4LogicalVolume* febCoolingLV = buildFEBCoolingBody(detectorGeo.getFEBCoolingGeometry());
1577 
1578  new G4PVPlacement(febCoolingTr, //Transformation
1579  febCoolingLV, //its logical volume
1580  "febCoolingLV", //its name
1581  detSupportLV, //its mother volume
1582  false, //no boolean operation
1583  iSlot); //copy number
1584  }
1585  }
1586 
1587  // G4LogicalVolume* hapdSupportPlateLV = new G4LogicalVolume(substraction, supportMaterial, "hapdSupport");
1588 
1589  G4Transform3D transform3 = G4Translate3D(0., 0., - backWallHeight / 2.);
1590  new G4PVPlacement(transform3, hapdSupportPlateLV, "supportPlate", detSupportLV, false, 1);
1591 
1592  // place electronics side neutron shield - for now hardcoded here, pending for more proper implementation!
1593  G4Box* shieldBox1 = new G4Box("shieldBox1", 20. / 2., 75. / 2., backWallHeight / 2.);
1594  G4Box* shieldBox2 = new G4Box("shieldBox2", 55. / 2., 40. / 2., backWallHeight / 2.);
1595  G4LogicalVolume* shield1 = new G4LogicalVolume(shieldBox1, Materials::get("BoratedPoly"), "ARICH.FWDShield1");
1596  G4LogicalVolume* shield2 = new G4LogicalVolume(shieldBox2, Materials::get("BoratedPoly"), "ARICH.FWDShield2");
1597  double dphi = 2 * M_PI / 36.;
1598  double r1 = wallR[0] - 15.;
1599  double r2 = wallR[0] - 15. - 20. / 2. - 55. / 2.;
1600  for (int i = 0; i < 36; i++) {
1601  double phi = (i + 0.5) * dphi;
1602  G4RotationMatrix rot;
1603  rot.rotateZ(phi);
1604  G4ThreeVector trans(r1 * cos(phi), r1 * sin(phi), detGeo.getSupportThickness() / 2.);
1605  G4ThreeVector trans1(r2 * cos(phi), r2 * sin(phi), detGeo.getSupportThickness() / 2.);
1606  new G4PVPlacement(G4Transform3D(rot, trans), shield1, "ARICH.FWDShield1", detSupportLV, false, i);
1607  new G4PVPlacement(G4Transform3D(rot, trans1), shield2, "ARICH.FWDShield2", detSupportLV, false, i);
1608  }
1609 
1610  return detSupportLV;
1611  }
1612 
1613 
1614  G4LogicalVolume* GeoARICHCreator::buildMirror(const ARICHGeometryConfig& detectorGeo)
1615  {
1616 
1617  const ARICHGeoMirrors& mirrGeo = detectorGeo.getMirrors();
1618 
1619  // read m_config
1620  string mirrMat = mirrGeo.getMaterial();
1621  G4Material* mirrorMaterial = Materials::get(mirrMat);
1622 
1623  double mThick = mirrGeo.getPlateThickness();
1624  double mLength = mirrGeo.getPlateLength();
1625  double mWidth = mirrGeo.getPlateWidth();
1626 
1627  G4Box* mirrPlate = new G4Box("mirrPlate", mThick / 2., mLength / 2., mWidth / 2.);
1628 
1629  G4LogicalVolume* lmirror = new G4LogicalVolume(mirrPlate, mirrorMaterial, "ARICH.mirrorPlate");
1630 
1631  Materials& materials = Materials::getInstance();
1632 
1633  G4OpticalSurface* optSurf = materials.createOpticalSurface(mirrGeo.getMirrorSurface());
1634  new G4LogicalSkinSurface("mirrorSurface", lmirror, optSurf);
1635 
1636  return lmirror;
1637 
1638  }
1639 
1640  double GeoARICHCreator::getAvgRINDEX(G4Material* material)
1641  {
1642  G4MaterialPropertiesTable* mTable = material->GetMaterialPropertiesTable();
1643  if (!mTable) return 0;
1644  G4MaterialPropertyVector* mVector = mTable->GetProperty("RINDEX");
1645  if (!mVector) return 0;
1646  G4bool b;
1647  return mVector->GetValue(2 * Unit::eV / Unit::MeV, b);
1648  }
1649 
1650  void GeoARICHCreator::makeJoint(G4Material* supportMaterial, const std::vector<double>& par, G4AssemblyVolume* assemblyWedge)
1651  {
1652 
1653  int size = par.size();
1654  if (size < 4 || size > 8) B2ERROR("GeoARICHCreator::makeJoint: invalid number of joint wedge parameters");
1655  double lenx = par.at(0);
1656  double leny = par.at(1);
1657  double lenz = par.at(2);
1658  double thick = par.at(3);
1659 
1660  G4Box* wedgeBox1 = new G4Box("wedgeBox1", thick / 2., lenx / 2., leny / 2.);
1661  G4Box* wedgeBox2 = new G4Box("wedgeBox2", lenz / 2., lenx / 2., thick / 2.);
1662 
1663  G4LogicalVolume* wedgeBox1LV = new G4LogicalVolume(wedgeBox1, supportMaterial, "ARICH.supportWedge");
1664  G4LogicalVolume* wedgeBox2LV = new G4LogicalVolume(wedgeBox2, supportMaterial, "ARICH.supportWedge");
1665 
1666  G4RotationMatrix Rm;
1667  G4ThreeVector Ta(0, 0, 0);
1668  G4Transform3D Tr;
1669  Tr = G4Transform3D(Rm, Ta);
1670 
1671  assemblyWedge->AddPlacedVolume(wedgeBox1LV, Tr);
1672 
1673  Ta.setX(lenz / 2. + thick / 2.);
1674  Ta.setZ(leny / 2. - thick / 2.);
1675  Tr = G4Transform3D(Rm, Ta);
1676  assemblyWedge->AddPlacedVolume(wedgeBox2LV, Tr);
1677 
1678  if (size == 4) return;
1679  double edge = par.at(4);
1680 
1681  G4Box* wedgeBox3 = new G4Box("wedgeBox3", lenz / 2., edge / 2., thick / 2.);
1682  G4LogicalVolume* wedgeBox3LV = new G4LogicalVolume(wedgeBox3, supportMaterial, "ARICH.supportWedge");
1683 
1684  Ta.setZ(leny / 2. - thick - thick / 2.);
1685  Tr = G4Transform3D(Rm, Ta);
1686  assemblyWedge->AddPlacedVolume(wedgeBox3LV, Tr);
1687 
1688  G4Trap* wedgeBoxTmp = new G4Trap("wedgeBoxTmp", thick, leny - 2 * thick, lenz, edge);
1689  G4LogicalVolume* wedgeBox4LV;
1690  if (size == 8) {
1691  G4Tubs* wedgeBoxTmp1 = new G4Tubs("wedgeBoxTmp1", 0.0, par.at(5), thick, 0, 2.*M_PI);
1692  G4RotationMatrix rotHole;
1693  G4ThreeVector transHole(par.at(6), par.at(7), 0);
1694  G4SubtractionSolid* wedgeBox4 = new G4SubtractionSolid("wedgeBox4", wedgeBoxTmp, wedgeBoxTmp1, G4Transform3D(rotHole, transHole));
1695  wedgeBox4LV = new G4LogicalVolume(wedgeBox4, supportMaterial, "ARICH.supportWedge");
1696  } else wedgeBox4LV = new G4LogicalVolume(wedgeBoxTmp, supportMaterial, "ARICH.supportWedge");
1697 
1698  Rm.rotateX(-M_PI / 2.);
1699  Ta.setX(thick / 2. + edge / 4. + lenz / 4.);
1700  Ta.setZ(-thick / 2. - edge / 2.);
1701  Tr = G4Transform3D(Rm, Ta);
1702  assemblyWedge->AddPlacedVolume(wedgeBox4LV, Tr);
1703 
1704  }
1705 
1706  // simple geometry for beamtest setup
1707  /* void GeoARICHCreator::createSimple(const GearDir& content, G4LogicalVolume& topVolume)
1708  {
1709 
1710  B2INFO("ARICH simple (beamtest) geometry will be built.")
1711  GearDir envelopeParams(content, "Envelope");
1712  double xSize = envelopeParams.getLength("xSize") / Unit::mm;
1713  double ySize = envelopeParams.getLength("ySize") / Unit::mm;
1714  double zSize = envelopeParams.getLength("zSize") / Unit::mm;
1715  string envMat = envelopeParams.getString("material");
1716  G4Material* envMaterial = Materials::get(envMat);
1717 
1718  G4Box* envBox = new G4Box("envBox", xSize / 2., ySize / 2., zSize / 2.);
1719  G4LogicalVolume* lenvBox = new G4LogicalVolume(envBox, envMaterial, "ARICH.envelope");
1720  new G4PVPlacement(G4Transform3D(), lenvBox, "ARICH.envelope", &topVolume, false, 1);
1721  setVisibility(*lenvBox, false);
1722  ARICHGeometryPar* m_arichgp = ARICHGeometryPar::Instance();
1723  m_arichgp->Initialize(content);
1724  GearDir moduleParam(content, "Detector/Module");
1725  G4LogicalVolume* detModule = buildModule(moduleParam);
1726 
1727  double detZpos = content.getLength("Detector/Plane/zPosition") / Unit::mm;
1728  double detThick = content.getLength("Detector/Module/moduleZSize") / Unit::mm;
1729  int nModules = m_arichgp->getNMCopies();
1730  for (int i = 1; i <= nModules; i++) {
1731  G4ThreeVector origin = m_arichgp->getOriginG4(i); origin.setZ(detZpos + detThick / 2.);
1732  double angle = m_arichgp->getModAngle(i);
1733  G4RotationMatrix Ra;
1734  Ra.rotateZ(angle);
1735  G4Transform3D trans = G4Transform3D(Ra, origin);
1736  new G4PVPlacement(G4Transform3D(Ra, origin), detModule, "detModule", lenvBox, false, i);
1737  }
1738 
1739  // place aerogel tiles
1740  GearDir aerogelParam(content, "Aerogel");
1741  double sizeX = aerogelParam.getLength("tileXSize") / Unit::mm;
1742  double sizeY = aerogelParam.getLength("tileYSize") / Unit::mm;
1743  double posX = aerogelParam.getLength("tileXPos") / Unit::mm;
1744  double posY = aerogelParam.getLength("tileYPos") / Unit::mm;
1745  int ilayer = 0;
1746  BOOST_FOREACH(const GearDir & layer, aerogelParam.getNodes("Layers/Layer")) {
1747  double posZ = layer.getLength("zPosition");
1748  double sizeZ = layer.getLength("thickness");
1749  string tileMat = layer.getString("material");
1750  G4Material* tileMaterial = Materials::get(tileMat);
1751  double rInd = getAvgRINDEX(tileMaterial);
1752  m_arichgp->setAeroRefIndex(ilayer, rInd);
1753  m_arichgp->setAerogelZPosition(ilayer, posZ);
1754  m_arichgp->setAerogelThickness(ilayer, sizeZ);
1755 
1756  G4MaterialPropertiesTable* mTable = tileMaterial->GetMaterialPropertiesTable();
1757  G4MaterialPropertyVector* mVector = mTable->GetProperty("RAYLEIGH");
1758  double lambda0 = 400;
1759  double e0 = 1240. / lambda0 * Unit::eV / Unit::MeV;
1760  G4bool b;
1761  m_arichgp->setAeroTransLength(ilayer, mVector->GetValue(e0, b)*Unit::mm);
1762  G4Box* tileBox = new G4Box("tileBox", sizeX / 2., sizeY / 2., sizeZ / 2. / Unit::mm);
1763  G4LogicalVolume* lTile = new G4LogicalVolume(tileBox, tileMaterial, "Tile", 0, m_sensitiveAero);
1764  setColor(*lTile, "rgb(0.0, 1.0, 1.0,1.0)");
1765  G4Transform3D trans = G4Translate3D(posX, posY, posZ / Unit::mm + sizeZ / 2. / Unit::mm);
1766  ilayer++;
1767  new G4PVPlacement(trans, lTile, "ARICH.tile", lenvBox, false, ilayer);
1768  }
1769 
1770  // place mirrors
1771  GearDir mirrorsParam(content, "Mirrors");
1772  double height = mirrorsParam.getLength("height") / Unit::mm;
1773  double width = mirrorsParam.getLength("width") / Unit::mm;
1774  double thickness = mirrorsParam.getLength("thickness") / Unit::mm;
1775  string mirrMat = mirrorsParam.getString("material");
1776  G4Material* mirrMaterial = Materials::get(mirrMat);
1777  G4Box* mirrBox = new G4Box("mirrBox", thickness / 2., height / 2., width / 2.);
1778  G4LogicalVolume* lmirror = new G4LogicalVolume(mirrBox, mirrMaterial, "mirror");
1779 
1780  Materials& materials = Materials::getInstance();
1781  GearDir surface(mirrorsParam, "Surface");
1782  G4OpticalSurface* optSurf = materials.createOpticalSurface(surface);
1783  new G4LogicalSkinSurface("mirrorsSurface", lmirror, optSurf);
1784  int iMirror = 0;
1785  BOOST_FOREACH(const GearDir & mirror, mirrorsParam.getNodes("Mirror")) {
1786  double xpos = mirror.getLength("xPos") / Unit::mm;
1787  double ypos = mirror.getLength("yPos") / Unit::mm;
1788  double zpos = mirror.getLength("zPos") / Unit::mm;
1789  double angle = mirror.getAngle("angle") / Unit::rad;
1790  G4ThreeVector origin(xpos, ypos, zpos + width / 2.);
1791  G4RotationMatrix Ra;
1792  Ra.rotateZ(angle);
1793  G4Transform3D trans = G4Transform3D(Ra, origin);
1794  new G4PVPlacement(G4Transform3D(Ra, origin), lmirror, "ARICH.mirror", lenvBox, false, iMirror);
1795  iMirror++;
1796  }
1797  }
1798  */
1799 
1800  } // namespace aricih
1802 } // namespace belleII
Geometry parameters of HAPD.
bool isSimple() const
Use simple aerogel configuration.
double getLayerThickness(unsigned iLayer) const
Get thickness of tiles i-th aerogel layer.
double getSupportThickness() const
Get support plate thickness.
unsigned getNRings() const
Get number of aluminum wall rings (should be number of tile rings + 1)
double getSupportOuterR() const
Get support plate outer radius.
double getCompensationARICHairVolumeThick_min() const
Get minimum thickness of the compensation volume with ARICH air.
double getImgTubeThickness() const
Get imaginary tube thikness just after aerogel layers used as volume to which tracks are extrapolated...
double getWallHeight() const
Get height of aluminum walls between aerogel tiles.
double getRingDPhi(unsigned iRing) const
Get phi (angle) distance between "phi" aluminum wall between aerogel tiles in i-th tile ring.
double getRotationY() const
Get angle of rotation around Y axis.
const std::string & getLayerMaterial(unsigned iLayer) const
Get material name of tiles i-th aerogel layer.
double getRotationZ() const
Get angle of rotation around Z axis.
const std::string & getSupportMaterial() const
Get material of support plate.
double getTileThickness(int ring, int column, int layerN) const
Get thickness of individual tile.
TVector3 getPosition() const
get position vector of aerogel plane in ARICH local frame
double getWallThickness() const
Get thickness of aluminum walls between aerogel tiles.
double getRingRadius(unsigned iRing) const
Get radius of i-th aluminum ring between aerogel tiles (inner radius of ring)
unsigned getNLayers() const
Get number of aerogel layers.
double getSupportInnerR() const
Get support plate inner radius.
int getFullAerogelMaterialDescriptionKey() const
Get full aerogel material description key.
std::string getTileMaterialName(int ring, int column, int layerN) const
Get material name of individual tile.
double getMaximumTotalTileThickness() const
Get maximum total thickness of the aerogel tiles tile_up + tile_down for all the slots.
double getTileGap() const
Get gap between aerogel tile and aluminum wall.
const std::vector< double > & getSimpleParams() const
Get parameters of simple aerogel configuration.
double getRotationX() const
Get angle of rotation around X axis.
Geometry parameters of cable envelope.
double getEnvelopeOuterRadius() const
Returns Outer radius of cables envelop.
const std::string & getCablesEffectiveMaterialName() const
Returns Effective material name describing cables.
TVector3 getEnvelopeCenterPosition() const
Returns position vector (TVector3) of cables envelop.
double getEnvelopeInnerRadius() const
Returns Inner radius of cables envelop.
double getEnvelopeThickness() const
Returns Thickness of cables envelop.
Geometry parameters of Cooling System.
double getColdTubeWallThickness() const
Returns cold tube wall thickness.
const std::vector< double > & getCoolingPosR() const
Returns vector of radial distance (r, pho) of the cooling system object center in polar coordinate sy...
const std::vector< double > & getCoolingL() const
Returns vector of lengs of the cooling system object with given geometry ID.
const std::string & getCoolingPipeMaterialName() const
Returns material name of cooling pipe.
double getEnvelopeOuterRadius() const
Returns Outer radius of cooling system assembly envelope.
double getColdTubeSpacing() const
Returns distance from center of the cold tube to edge of cooling plate.
const std::vector< double > & getCoolingPosPhi() const
Returns vector of azimuthal angle of the cooling system object center in polar coordinate system in d...
double getColdTubeR() const
Returns radius of cold tubes.
const std::vector< double > & getCoolingGeometryID() const
Returns vector of cooling system object geometry ID.
const std::string & getCoolingTestPlateMaterialName() const
Returns material name of cooling test plates.
TVector3 getEnvelopeCenterPosition() const
Returns position vector (TVector3) of cooling system assembly envelope.
TVector3 getCoolingTestPlateslengths() const
Returns sizes vector (TVector3) of cooling test plates.
double getRmax() const
Returns Size of cooling system pipe : outer radius in mm.
const std::string & getColdTubeMaterialName() const
Returns material name of cold tube.
const std::vector< double > & getCoolinRotationAngle() const
Returns vector of azimuthal angle of rotation aroud Z - axis of the cooling system object in polar co...
double getEnvelopeInnerRadius() const
Returns Inner radius of cooling system assembly envelope.
double getRmin() const
Returns Size of cooling system pipe : inner radius in mm.
double getDepthColdTubeInPlate() const
Returns depth of the cold tube in the cooling plate.
double getEnvelopeThickness() const
Returns Thickness of cooling system assembly envelope.
double getColdTubeSubtractedR() const
Returns outer radius of subtracted tubes for cold tube.
int getColdTubeNumber() const
Returns number of cold tubes in one plate.
Geometry parameters of ARICH photon detector plane.
unsigned getNSlots() const
Get total number of module slots.
unsigned getSlotRing(unsigned slotID) const
Get ring number of slot with slot ID.
double getSupportBackWallHeight() const
Get height of the aluminum walls between modules on the electronics side of aluminum support plate.
double getSlotR(unsigned modID) const
Get radial position of module with given module ID number.
double getSupportBackWallThickness() const
Get thickness of aluminum walls between modules on the electronics side of aluminum support plate.
double getSupportThickness() const
Get support plate thickness.
unsigned getNRings() const
Get number of module slot rings.
double getSupportOuterR() const
Get support plate outer radius.
double getSupportZPosition() const
Get Z position of support plate (start point in Z)
double getRingDPhi(unsigned iRing) const
Get phi (angle) distance between module slots in i-th ring.
double getRotationY() const
Get angle of rotation around Y axis.
double getRotationZ() const
Get angle of rotation around Z axis.
double getSlotPhi(unsigned modID) const
Get phi (angle) position of module with given module ID number.
double getModuleHoleSize() const
Get size of module hole in support plate.
const std::string & getSupportMaterial() const
Get material of support plate.
TVector3 getPosition() const
Get center point.
double getSupportInnerR() const
Get support plate inner radius.
double getRingR(unsigned iRing) const
Get radius of i-th module slot ring (radius of center point)
double getRotationX() const
Get angle of rotation around X axis.
Geometry parameters of Cooling System - version2 (v2).
double getBigSquareThickness() const
Returns thickness of the big square in mm.
double getBigSquareSize() const
Returns size of the big square in mm.
double getSmallSquareSize() const
Returns size of the small square in mm.
double getRectangleThickness() const
Returns thickness of the rectangle in mm.
const std::vector< double > & getFebcoolingv2GeometryID() const
Returns vector of feb cooling configuration/geometry ID.
double getSmallSquareThickness() const
Returns thickness of the small square in mm.
double getRectangleW() const
Returns width of the rectangle in mm.
double getRectangleL() const
Returns length of the rectangle in mm.
double getRectangleDistanceFromCenter() const
Returns distance from center of the rectangle in mm.
double getGamma() const
Returns rotation angle around z.
double getX() const
Returns translation in x.
double getAlpha() const
Returns rotation angle around x.
double getZ() const
Returns translation in z.
double getY() const
Returns translation in y.
double getBeta() const
Returns rotation angle around y.
Geometry parameters of HAPD.
Definition: ARICHGeoHAPD.h:24
double getSizeZ() const
Returns HAPD size in z.
Definition: ARICHGeoHAPD.h:174
const std::string & getFillMaterial() const
Returns fill (inside) material name.
Definition: ARICHGeoHAPD.h:271
const GeoOpticalSurface & getAPDSurface() const
Returns APD reflective optical surface.
Definition: ARICHGeoHAPD.h:295
const std::string & getWinMaterial() const
Returns window material name.
Definition: ARICHGeoHAPD.h:277
double getAPDSizeZ() const
Returns APD size in z.
Definition: ARICHGeoHAPD.h:204
double getWinThickness() const
Returns window thickness.
Definition: ARICHGeoHAPD.h:186
const std::string & getWallMaterial() const
Returns wall (casing) material name.
Definition: ARICHGeoHAPD.h:265
const std::string & getAPDMaterial() const
Returns APD material name.
Definition: ARICHGeoHAPD.h:283
const std::string & getFEBMaterial() const
Returns FEB material name.
Definition: ARICHGeoHAPD.h:289
double getSizeX() const
Returns HAPD size in x.
Definition: ARICHGeoHAPD.h:162
double getSizeY() const
Returns HAPD size in y.
Definition: ARICHGeoHAPD.h:168
double getFEBSizeY() const
Returns FEB size in y.
Definition: ARICHGeoHAPD.h:216
double getWallThickness() const
Returns wall thickness.
Definition: ARICHGeoHAPD.h:180
double getFEBSizeZ() const
Returns FEB size in z.
Definition: ARICHGeoHAPD.h:222
double getModuleSizeZ() const
Returns module size in z (HAPD + FEB height)
Definition: ARICHGeoHAPD.h:228
double getAPDSizeY() const
Returns APD size in y.
Definition: ARICHGeoHAPD.h:198
double getFEBSizeX() const
Returns FEB size in x.
Definition: ARICHGeoHAPD.h:210
double getAPDSizeX() const
Returns APD size in x.
Definition: ARICHGeoHAPD.h:192
double getRotationY() const
Get angle of rotation around Y axis.
double getRotationZ() const
Get angle of rotation around Z axis.
double getOuterRadius() const
Get ARICH master volume outer radius.
double getInnerRadius() const
Get ARICH master volume inner radius.
TVector3 getPosition() const
Get position of ARICH master volume center point in global Belle II coordinates.
const std::string & getMaterial() const
Get material of ARICH master volume.
double getLength() const
Get ARICH master volume length.
double getRotationX() const
Get angle of rotation around X axis.
Geometry parameters of Merger PCB.
const std::vector< double > & getMergerSlotID() const
Returns vector of merger boards slot numbers.
double getMergerPCBscrewholeR() const
Returns merger PCB screw hole radius.
double getSingleMergerEnvelopeThickness() const
Returns single merger PCB and merger cooling envelop thickness.
double getSizeW() const
Returns merger PCB width.
const std::string & getMergerPCBMaterialName() const
Returns merger PCB material name.
double getMergerPCBscrewholePosdX1() const
Returns merger PCB screw hole position from the bottom edge.
double getEnvelopeOuterRadius() const
Returns Outer radius of merger PCB assembly envelope.
double getMergerPCBscrewholePosdX2() const
Returns merger PCB screw hole position from the top edge.
TVector3 getEnvelopeCenterPosition() const
Returns position vector (TVector3) of merger PCB assembly envelope.
TVector3 getSingleMergeEnvelopePosition() const
Returns position vector (TVector3) of merger PCB inside the single merger envelope.
double getThickness() const
Returns merger PCB thickness.
double getSingleMergerEnvelopeSizeW() const
Returns single merger PCB and merger cooling envelop width.
double getEnvelopeInnerRadius() const
Returns Inner radius of merger PCB assembly envelope.
const std::vector< double > & getMergerPosR() const
Returns vector of merger boards distances from the center in mm.
const std::vector< double > & getMergerOrientation() const
Returns vector of merger boarts orientations in deg.
double getSingleMergerEnvelopeSizeL() const
Returns single merger PCB and merger cooling envelop length.
double getSizeL() const
Returns merger PCB lenght.
double getEnvelopeThickness() const
Returns Thickness of merger PCB assembly envelope.
double getMergerPCBscrewholePosdY() const
Returns merger PCB screw hole position from the left and right sides.
const std::vector< double > & getSingleMergerenvelopeDeltaZ() const
Returns vector of Z position of the single merger and merger cooling body envelope inside global merg...
const std::vector< double > & getMergerAngle() const
Returns vector of merger boarts azimuthal angles in polar coordinate system in deg.
const ARICHPositionElement & getDisplacementElement(int mirrorID) const
Returns displacement parameters for given mirror plate.
Geometry parameters of HAPD.
double getZPosition() const
Get nominal Z position of mirror plates (center point in ARICH local frame)
const GeoOpticalSurface & getMirrorSurface() const
Returns mirror reflective optical surface.
unsigned getNMirrors() const
Get number of mirror plates.
double getPlateWidth() const
Get width of mirror plate.
const std::string & getMaterial() const
Get material name of mirror plates.
double getPlateLength() const
Get length of mirror plate.
double getRadius() const
Get nominal radius at which mirror plates are placed (center of plate)
double getStartAngle() const
Get phi angle of position of the first mirror plate.
double getPlateThickness() const
Get thickness of mirror plate.
Geometry parameters of ARICH support structures and neutron shield.
double getWedgeR(unsigned i) const
Get radius at which i-th wedge is placed.
double getWedgePhi(unsigned i) const
Get phi angle at which i-th wedge is placed.
const std::string & getWedgeMaterial(unsigned i) const
Get material of i-th wedge.
const std::string & getTubeMaterial(unsigned i) const
Get material of i-th tube.
double getTubeInnerR(unsigned i) const
Get tube inner radius.
double getTubeLength(unsigned i) const
Get tube length.
double getTubeOuterR(unsigned i) const
Get tube outer radius.
double getWedgeZ(unsigned i) const
Get Z position of i-th wedge.
unsigned getNTubes() const
Get number of tube volumes to be placed.
const std::vector< double > getWedge(unsigned i) const
Get parameters of wedge.
int getWedgeType(unsigned i) const
Get type of i-th wedge.
const std::string & getTubeName(unsigned i) const
Get name of i-th tube.
unsigned getNWedges() const
Get number of wedges to be placed.
double getTubeZPosition(unsigned i) const
Get tube Z position.
box getBox(unsigned i) const
Get box paramaters.
unsigned getNBoxes() const
Get number of box volumes.
The Class for ARICH Geometry Parameters.
const ARICHGeoMirrorDisplacement & getMirrorDisplacement() const
Get mirror displacement parameters.
const ARICHGeoCooling & getCoolingGeometry() const
Get ARICH cooling system geometry parameters.
const ARICHGeoMirrors & getMirrors() const
Get mirrors geometry configuration.
const ARICHGeoCablesEnvelope & getCablesEnvelope() const
Get ARICH cables envelop geometry parameters.
static void useGeantUnits()
Use Geant4 units when returning geometry parameters.
const ARICHGeoHAPD & getHAPDGeometry() const
Get HAPD geometry parameters.
const ARICHGeoGlobalDisplacement & getGlobalDisplacement() const
Get global displacement parameters.
const ARICHGeoFEBCooling & getFEBCoolingGeometry() const
Get ARICH FEB cooling system (v2) geometry parameters.
const ARICHGeoMerger & getMergerGeometry() const
Get Merger PCB geometry parameters.
const ARICHGeoAerogelPlane & getAerogelPlane() const
Get geometry configuration of aerogel plane.
const ARICHGeoSupport & getSupportStructure() const
Get ARICH support structure geometry configuration.
static void useBasf2Units()
Use basf2 units when returning geometry parameters.
const ARICHGeoDetectorPlane & getDetectorPlane() const
Get geometry configuration of HAPD plane.
const ARICHGeoMasterVolume & getMasterVolume() const
Get ARICH master volume geometry configuration.
Position element for ARICH.
The Class for BeamBackground Sensitive Detector.
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
G4LogicalVolume * buildSimpleAerogelPlane(const ARICHGeometryConfig &detectorGeo)
build simple aerogel plane (for cosmic test)
G4LogicalVolume * buildHAPD(const ARICHGeoHAPD &hapdPar)
build the HAPD modules
void makeJoint(G4Material *supportMaterial, const std::vector< double > &pars, G4AssemblyVolume *assemblyWedge)
build joints of the ARICH support structure
virtual ~GeoARICHCreator()
The destructor of the GeoARICHreator class.
OptionalDBObjPtr< ARICHGeoMergerCooling > m_mergerCooling
merger cooling bodies geometry from the DB
G4LogicalVolume * buildMergerPCBEnvelopePlane(const ARICHGeometryConfig &detectorGeo)
build merger PCB assembly envelope plane
G4LogicalVolume * buildFEBCoolingBody(const ARICHGeoFEBCooling &coolingv2Geo)
build FEB cooling bodies (cooling system update after phase 2)
G4LogicalVolume * buildCoolingTorus(const unsigned i_volumeID, const ARICHGeoCooling &coolingGeo)
build cooling tube (G4Torus)
G4LogicalVolume * buildAerogelPlaneAveragedOverLayers(const ARICHGeometryConfig &detectorGeo)
build aerogel plane with average properties of aerogel per layer
G4LogicalVolume * buildCables(const ARICHGeoCablesEnvelope &cablesGeo)
build the cables envelop with effective material describing cables
G4LogicalVolume * buildMergerCooling(unsigned iType)
build merger cooling bodies (cooling system update after phase 2)
double getAvgRINDEX(G4Material *material)
get refractive index of the material
G4LogicalVolume * buildCoolingTube(const unsigned i_volumeID, const ARICHGeoCooling &coolingGeo)
build cooling tube (G4Tubs)
void createGeometry(G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Create detector geometry.
G4LogicalVolume * buildDetectorSupportPlate(const ARICHGeometryConfig &detectorGeo)
build detector support plate
DBObjPtr< ARICHModulesInfo > m_modInfo
information on installed modules from the DB
G4LogicalVolume * buildAerogelPlane(const ARICHGeometryConfig &detectorGeo)
build aerogel plane
SensitiveAero * m_sensitiveAero
pointer to sensitive aerogel - used instead of tracking
G4LogicalVolume * buildDetectorPlane(const ARICHGeometryConfig &detectorGeo)
build detector plane
ARICHGeometryConfig m_config
geometry configuration
G4LogicalVolume * buildCoolingEnvelopePlane(const ARICHGeoCooling &coolingGeo)
build cooling system assembly envelope plane
int m_isBeamBkgStudy
flag the beam background study
G4LogicalVolume * buildMergerEnvelope(const ARICHGeoMerger &mergerGeo, int type)
build single merger and merger cooling body envelope logical volume
G4LogicalVolume * buildMerger(const ARICHGeoMerger &mergerGeo)
build the merger PCB logical volume
SensitiveDetector * m_sensitive
pointer to sensitive detector
G4LogicalVolume * buildAerogelPlaneWithIndividualTilesProp(const ARICHGeometryConfig &detectorGeo)
with individual properties of aerogel tiles
G4LogicalVolume * buildMirror(const ARICHGeometryConfig &detectorGeo)
build mirrors
G4LogicalVolume * buildCoolingTestPlate(const ARICHGeoCooling &coolingGeo)
build cooling test plates
This is optional (temporary) class that provides information on track parameters on aerogel plane,...
Definition: SensitiveAero.h:22
The Class for ARICH Sensitive Detector.
Thin wrapper around the Geant4 Material system.
Definition: Materials.h:48
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
static Materials & getInstance()
Get a reference to the singleton instance.
Definition: Materials.cc:85
G4OpticalSurface * createOpticalSurface(const gearbox::Interface &parameters)
Create an optical surface from parameters, will abort on error.
Definition: Materials.cc:291
int doBeamBackgroundStudy() const
returns 1 if beam background study (to add additional sensitive modules, detect neutrons,...
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
void setVisibility(G4LogicalVolume &volume, bool visible)
Helper function to quickly set the visibility of a given volume.
Definition: utilities.cc:108
void setColor(G4LogicalVolume &volume, const std::string &color)
Set the color of a logical volume.
Definition: utilities.cc:100
GeometryTypes
Flag indiciating the type of geometry to be used.
Abstract base class for different kinds of events.
Struct to hold parameters of box volumes (examples, scintilators for cosmic test)
Structure which holds apexes of the tessellation volumes.