Belle II Software  release-06-00-14
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 = makeJoint(wedge1Material, wedge1);
315  G4AssemblyVolume* assemblyWedge2 = makeJoint(wedge2Material, wedge2);
316 
317  G4Transform3D Tr;
318  int nWedge = supportPar.getNWedges();
319  for (int i = 0; i < nWedge; i++) {
320  G4RotationMatrix Rx;
321  int wgtype = supportPar.getWedgeType(i);
322  double phi = supportPar.getWedgePhi(i);
323  Rx.rotateZ(phi);
324  double r = supportPar.getWedgeR(i);
325  double z = supportPar.getWedgeZ(i);
326  G4ThreeVector Tb(r * cos(phi), r * sin(phi), z);
327 
328  if (shieldLV[0]) {
329  z -= shieldZ[0];
330  Tb.setZ(z);
331  Tr = G4Transform3D(Rx, Tb);
332  if (wgtype == 1) assemblyWedge1->MakeImprint(shieldLV[0], Tr);
333  else if (wgtype == 2) assemblyWedge2->MakeImprint(shieldLV[0], Tr);
334  continue;
335  }
336 
337  Tr = G4Transform3D(Rx, Tb);
338  if (wgtype == 1) assemblyWedge1->MakeImprint(masterLV, Tr);
339  else if (wgtype == 2) assemblyWedge2->MakeImprint(masterLV, Tr);
340  else B2ERROR("GeoARICHCreator: invalid support wedge type!");
341  }
342 
343  // place neutron shield volumes
344  if (shieldLV[0])
345  new G4PVPlacement(G4Transform3D(G4RotationMatrix(), G4ThreeVector(0, 0, shieldZ[0])), shieldLV[0], "ARICH.NeutronShield1", masterLV,
346  false, 1);
347  if (shieldLV[1])
348  new G4PVPlacement(G4Transform3D(G4RotationMatrix(), G4ThreeVector(0, 0, shieldZ[1])), shieldLV[1], "ARICH.NeutronShield2", masterLV,
349  false, 1);
350 
351  // place mirror holders, for now skip
352  /*
353  double mirSupX = 17.;
354  double mirSupY = 14.;
355  double mirSupLen = 131;
356  double mirSupLen1 = 201;
357 
358  std::vector<G4TwoVector> polygon;
359  polygon.assign(12, G4TwoVector());
360  polygon[11] = G4TwoVector(-mirSupX / 2., -mirSupY / 2.);
361  polygon[10] = G4TwoVector(-mirSupX / 2., 2.);
362  polygon[9] = G4TwoVector(-mirSupX / 4., 2. - 0.75);
363  polygon[8] = G4TwoVector(-mirSupX / 4., 5. - 0.75);
364  polygon[7] = G4TwoVector(-mirSupX / 2., 5);
365  polygon[6] = G4TwoVector(-mirSupX / 2., mirSupY / 2.);
366  polygon[5] = G4TwoVector(mirSupX / 2., mirSupY / 2.);
367  polygon[4] = G4TwoVector(mirSupX / 2., 6);
368  polygon[3] = G4TwoVector(mirSupX / 4., 6. - 0.75);
369  polygon[2] = G4TwoVector(mirSupX / 4., 2. - 0.75);
370  polygon[1] = G4TwoVector(mirSupX / 2., 2.);
371  polygon[0] = G4TwoVector(mirSupX / 2., -mirSupY / 2.);
372 
373  std::vector<G4ExtrudedSolid::ZSection> zsections;
374  zsections.push_back(G4ExtrudedSolid::ZSection(-mirSupLen / 2., G4TwoVector(0, 0), 1));
375  zsections.push_back(G4ExtrudedSolid::ZSection(mirSupLen / 2., G4TwoVector(0, 0), 1));
376 
377  G4ExtrudedSolid* shape = new G4ExtrudedSolid("bla", polygon, zsections);
378  G4LogicalVolume* shapeLV = new G4LogicalVolume(shape, wedge1Material, "mirrorsupport");
379 
380  G4Box* shape1 = new G4Box("shape1", mirSupX / 2., 2. / 2., mirSupLen1 / 2.);
381  G4LogicalVolume* shape1LV = new G4LogicalVolume(shape1, wedge1Material, "mirrorholder");
382 
383  G4AssemblyVolume* mirroHolder = new G4AssemblyVolume();
384  G4RotationMatrix Rh;
385  G4ThreeVector Th(0, 0, 0);
386  Tr = G4Transform3D(Rh, Th);
387  mirroHolder->AddPlacedVolume(shapeLV, Tr);
388 
389  Th.setZ(-35.0);
390  Th.setY(-mirSupY / 2. - 1.);
391  Tr = G4Transform3D(Rh, Th);
392  mirroHolder->AddPlacedVolume(shape1LV, Tr);
393 
394  double rholder = m_config.getDetectorPlane().getSupportOuterR() - mirSupY / 2.;
395 
396  angl = dphi / 2.;
397  for (int i = 1; i < nMirrors + 1; i++) {
398  G4RotationMatrix rotMirror;
399  rotMirror.rotateX(M_PI);
400  rotMirror.rotateZ(-M_PI / 2. + angl);
401  G4ThreeVector transMirror(rholder * cos(angl), rholder * sin(angl), mirZPos + zShift);
402  Tr = G4Transform3D(rotMirror, transMirror);
403  mirroHolder->MakeImprint(masterLV, Tr);
404  angl += dphi;
405  }
406  */
407 
408  // temporary for simple cosmic test
409  // if using simple configuration, place scinitilators
411  int nBoxes = supportPar.getNBoxes();
412  for (int i = 0; i < nBoxes; i++) {
413  ARICHGeoSupport::box box = supportPar.getBox(i);
414  G4Box* scintBox = new G4Box("scintBox", box.size[0] * 10. / 2., box.size[1] * 10. / 2., box.size[2] * 10. / 2.);
415  G4Material* scintMaterial = Materials::get(box.material);
416  G4LogicalVolume* scintLV = new G4LogicalVolume(scintBox, scintMaterial, box.name);
417  scintLV->SetSensitiveDetector(m_sensitiveAero);
418  G4RotationMatrix rotScint;
419  rotScint.rotateX(box.rotation[0]);
420  rotScint.rotateY(box.rotation[1]);
421  rotScint.rotateZ(box.rotation[2]);
422  TVector3 transScintTV(box.position[0], box.position[1], box.position[2]);
423  transScintTV = m_config.getMasterVolume().pointToGlobal(transScintTV);
424  B2INFO("GeoARICHCreator: Scintilator " << box.name << " placed at global: " << transScintTV.X() << " " << transScintTV.Y() << " " <<
425  transScintTV.Z());
426  G4ThreeVector transScint(box.position[0] * 10., box.position[1] * 10., box.position[2] * 10.);
427  new G4PVPlacement(G4Transform3D(rotScint, transScint), scintLV, "scintilator", masterLV, false, 1);
428  }
429  }
430 
432 
433  return;
434 
435  }
436 
437 
439  {
440 
441  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
442 
443  const std::vector<double>& params = aeroGeo.getSimpleParams();
444 
445  // support plane
446  double rin = aeroGeo.getSupportInnerR();
447  double rout = aeroGeo.getSupportOuterR();
448  double thick = aeroGeo.getSupportThickness();
449  string supportMat = aeroGeo.getSupportMaterial();
450  double wallHeight = aeroGeo.getWallHeight();
451  G4Material* supportMaterial = Materials::get(supportMat);
452  G4Material* gapMaterial = Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
453 
454  // master volume
455  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight) / 2., 0, 2 * M_PI);
456  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
457 
458  // support plate
459  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
460  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
461  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
462 
463  unsigned nLayer = aeroGeo.getNLayers();
464 
465  // loop over layers
466  double zLayer = 0;
467  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
468  double layerThick = aeroGeo.getLayerThickness(iLayer);
469  std::stringstream tileName;
470  tileName << "aerogelTile_" << iLayer;
471 
472  G4Box* tileShape = new G4Box(tileName.str(), params[0] * 10. / 2., params[1] * 10. / 2., layerThick / 2.);
473 
474  G4Material* aeroMaterial = Materials::get(aeroGeo.getLayerMaterial(iLayer));
475 
476  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, aeroMaterial, string("ARICH.") + tileName.str());
477 
478  G4ThreeVector transTile(params[2] * 10., params[3] * 10., (thick + layerThick - wallHeight) / 2. + zLayer);
479  G4RotationMatrix Ra;
480  Ra.rotateZ(params[4]);
481  new G4PVPlacement(G4Transform3D(Ra, transTile), tileLV, string("ARICH.") + tileName.str(), aerogelPlaneLV, false, iLayer);
482 
483  zLayer += layerThick;
484  }
485 
486  new G4PVPlacement(G4Translate3D(0., 0., -wallHeight / 2.), supportTubeLV, "ARICH.AerogelSupportPlate", aerogelPlaneLV, false, 1);
487 
488  return aerogelPlaneLV;
489  }
490 
491  G4LogicalVolume* GeoARICHCreator::buildAerogelPlane(const ARICHGeometryConfig& detectorGeo)
492  {
493  if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 0)
494  return buildAerogelPlaneAveragedOverLayers(detectorGeo);
495  else if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 1)
496  return buildAerogelPlaneWithIndividualTilesProp(detectorGeo);
497  else
498  B2ERROR("GeoARICHCreator::buildAerogelPlane --> getFullAerogelMaterialDescriptionKey() is wrong");
499  return NULL;
500  }
501 
503  {
504 
505  //cout<<"GeoARICHCreator::buildAerogelPlaneAveragedOverLayers(const ARICHGeometryConfig& detectorGeo)"<<endl;
506 
507  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
508 
509  // support plane
510  double rin = aeroGeo.getSupportInnerR();
511  double rout = aeroGeo.getSupportOuterR();
512  double thick = aeroGeo.getSupportThickness();
513  double wallThick = aeroGeo.getWallThickness();
514  double wallHeight = aeroGeo.getWallHeight();
515  string supportMat = aeroGeo.getSupportMaterial();
516  G4Material* supportMaterial = Materials::get(supportMat);
517  G4Material* gapMaterial = Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
518  G4Material* imgMaterial = Materials::get("ARICH_Air");
519  // master volume
520 
521  //double imgTubeLen = 0.5; // if changed, change position of aerogel plane also in main function!!
522  double imgTubeLen = aeroGeo.getImgTubeThickness();
523  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight + imgTubeLen) / 2., 0, 2 * M_PI);
524  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
525 
526  // support plate
527  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
528  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
529  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
530 
531  // imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
532  G4Tubs* imgTube = new G4Tubs("imgTube", rin, rout, imgTubeLen / 2., 0, 2 * M_PI);
533  G4LogicalVolume* imgTubeLV = new G4LogicalVolume(imgTube, imgMaterial, "ARICH.AerogelImgPlate");
534  imgTubeLV->SetSensitiveDetector(m_sensitiveAero);
535 
536 
537  // read radiuses of aerogel slot aluminum walls
538  std::vector<double> wallR;
539  unsigned nRing = aeroGeo.getNRings();
540  for (unsigned iRing = 1; iRing < nRing + 1; iRing++) {
541  wallR.push_back(aeroGeo.getRingRadius(iRing));
542  }
543 
544  unsigned nLayer = aeroGeo.getNLayers();
545  double tileGap = aeroGeo.getTileGap();
546  G4Transform3D transform = G4Translate3D(0., 0., (thick - imgTubeLen) / 2.);
547 
548  for (unsigned iRing = 0; iRing < nRing; iRing++) {
549 
550  // aluminum walls between tile rings (r wall)
551  std::stringstream wallName;
552  wallName << "supportWallR_" << iRing + 1;
553  G4Tubs* supportWall = new G4Tubs(wallName.str().c_str(), wallR[iRing], wallR[iRing] + wallThick, wallHeight / 2., 0, 2 * M_PI);
554  G4LogicalVolume* supportWallLV = new G4LogicalVolume(supportWall, supportMaterial, string("ARICH.") + wallName.str().c_str());
555 
556  new G4PVPlacement(transform, supportWallLV, string("ARICH.") + wallName.str().c_str(), aerogelPlaneLV, false, 0);
557 
558  if (iRing == 0) continue;
559 
560  // place phi aluminum walls
561  double dphi = aeroGeo.getRingDPhi(iRing);
562 
563  wallName.str("");
564  wallName << "supportWallPhi_" << iRing + 1;
565  G4Box* wall = new G4Box(wallName.str(), (wallR[iRing] - wallR[iRing - 1] - wallThick) / 2. - 1., thick / 2., wallHeight / 2.);
566  G4LogicalVolume* wallLV = new G4LogicalVolume(wall, supportMaterial, string("ARICH.") + wallName.str());
567  double r = (wallR[iRing - 1] + wallThick + wallR[iRing]) / 2.;
568  double zLayer = 0;
569 
570  // loop over layers
571  int iSlot = 1;
572  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
573  double iphi = 0;
574  double layerThick = aeroGeo.getLayerThickness(iLayer);
575 
576  std::stringstream tileName;
577  tileName << "aerogelTile_" << iRing << "_" << iLayer;
578 
579  G4Tubs* tileShape = new G4Tubs(tileName.str(), wallR[iRing - 1] + wallThick + tileGap, wallR[iRing] - tileGap, layerThick / 2.,
580  (tileGap + wallThick / 2.) / wallR[iRing], dphi - (2.*tileGap + wallThick) / wallR[iRing]);
581 
582  G4Material* aeroMaterial = Materials::get(aeroGeo.getLayerMaterial(iLayer));
583  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, aeroMaterial, string("ARICH.") + tileName.str());
584 
585  while (iphi < 2 * M_PI - 0.0001) {
586  G4ThreeVector trans(r * cos(iphi), r * sin(iphi), (thick - imgTubeLen) / 2.);
587  G4RotationMatrix Ra;
588  Ra.rotateZ(iphi);
589 
590  if (iLayer == 1) new G4PVPlacement(G4Transform3D(Ra, trans), wallLV, string("ARICH.") + wallName.str(), aerogelPlaneLV, false,
591  iSlot);
592 
593  G4ThreeVector transTile(0, 0, (thick + layerThick - wallHeight - imgTubeLen) / 2. + zLayer);
594  new G4PVPlacement(G4Transform3D(Ra, transTile), tileLV, string("ARICH.") + tileName.str(), aerogelPlaneLV, false, iSlot);
595  iphi += dphi;
596  iSlot++;
597  }
598  zLayer += layerThick;
599  }
600  }
601 
602  new G4PVPlacement(G4Translate3D(0., 0., -(wallHeight + imgTubeLen) / 2.), supportTubeLV, "ARICH.AerogelSupportPlate",
603  aerogelPlaneLV,
604  false, 1);
605 
606  new G4PVPlacement(G4Translate3D(0., 0., (wallHeight + thick) / 2.), imgTubeLV, "ARICH.AerogelImgPlate", aerogelPlaneLV, false, 1);
607 
608  return aerogelPlaneLV;
609 
610  }
611 
613  {
614 
615  //cout << "GeoARICHCreator::buildAerogelPlaneWithIndividualTilesProp(const ARICHGeometryConfig& detectorGeo)" << endl;
616 
617  const ARICHGeoAerogelPlane& aeroGeo = detectorGeo.getAerogelPlane();
618 
619  // support plane
620  double rin = aeroGeo.getSupportInnerR();
621  double rout = aeroGeo.getSupportOuterR();
622  double thick = aeroGeo.getSupportThickness();
623  double wallThick = aeroGeo.getWallThickness();
624  // Maximum total (up and down) thickness of the aerogel tiles
625  double maxTotalAerogelThick = aeroGeo.getMaximumTotalTileThickness();
626  //cout<<"maxTotalAerogelThick "<<maxTotalAerogelThick<<endl
627  // <<"wallHeight "<<wallHeight<<endl;
628  // In case of individual thickness of the tiles we need to define compensation
629  // volume with ARICH air. This volume situated between aerogel tile and image plane (imgTube).
630  // Minimum thickness of the compensation volume with ARICH air
631  double compensationARICHairVolumeThick_min = aeroGeo.getCompensationARICHairVolumeThick_min(); // mm
632  // Please note redefinition of wallHeight value
633  double wallHeight = maxTotalAerogelThick + compensationARICHairVolumeThick_min;
634 
635  string supportMat = aeroGeo.getSupportMaterial();
636  G4Material* supportMaterial = Materials::get(supportMat);
637 
638  G4Material* gapMaterial =
639  Materials::get("Air"); // Air without refractive index (to kill photons, to mimic black paper around tile)
640  G4Material* imgMaterial = Materials::get("ARICH_Air"); // Air with defined optical properties to propagate cherenkov photons
641 
642  // master volume
643  double imgTubeLen = aeroGeo.getImgTubeThickness(); // if changed, change position of aerogel plane also in main function!!
644  G4Tubs* aerogelTube = new G4Tubs("aerogelTube", rin, rout, (thick + wallHeight + imgTubeLen) / 2., 0, 2 * M_PI);
645  G4LogicalVolume* aerogelPlaneLV = new G4LogicalVolume(aerogelTube, gapMaterial, "ARICH.AaerogelPlane");
646 
647  // support plate
648  G4Tubs* supportTube = new G4Tubs("aeroSupportTube", rin, rout, thick / 2., 0, 2 * M_PI);
649  G4LogicalVolume* supportTubeLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.AerogelSupportPlate");
650  //supportTubeLV->SetSensitiveDetector(m_sensitiveAero);
651 
652  // imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
653  G4Tubs* imgTube = new G4Tubs("imgTube", rin, rout, imgTubeLen / 2., 0, 2 * M_PI);
654  G4LogicalVolume* imgTubeLV = new G4LogicalVolume(imgTube, imgMaterial, "ARICH.AerogelImgPlate");
655  imgTubeLV->SetSensitiveDetector(m_sensitiveAero);
656 
657  // read radiuses of aerogel slot aluminum walls
658  std::vector<double> wallR;
659  unsigned nRing = aeroGeo.getNRings();
660  for (unsigned iRing = 1; iRing < nRing + 1; iRing++) {
661  wallR.push_back(aeroGeo.getRingRadius(iRing));
662  }
663 
664  unsigned nLayer = aeroGeo.getNLayers();
665  double tileGap = aeroGeo.getTileGap();
666  G4Transform3D transform = G4Translate3D(0., 0., (thick - imgTubeLen) / 2.);
667 
668  for (unsigned iRing = 0; iRing < nRing; iRing++) {
669 
670  // Aluminum walls between tile rings (r wall)
671  std::stringstream wallName;
672  wallName << "supportWallR_" << iRing + 1;
673  //cout<<"wallName = "<<wallName.str().c_str()<<endl;
674  G4Tubs* supportWall = new G4Tubs(wallName.str().c_str(), wallR[iRing], wallR[iRing] + wallThick, wallHeight / 2., 0, 2 * M_PI);
675  G4LogicalVolume* supportWallLV = new G4LogicalVolume(supportWall, supportMaterial, string("ARICH.") + wallName.str().c_str());
676  new G4PVPlacement(transform, supportWallLV, string("ARICH.") + wallName.str().c_str(), aerogelPlaneLV, false, 0);
678 
679  // There are only 4 rings of aerogel - the first one is only for mechanical support
680  if (iRing == 0) continue;
681 
682  // dphi - distance between centers of two consecutive aluminum walls (diaphragm) in one ring
683  double dphi = aeroGeo.getRingDPhi(iRing);
684 
685  // Aluminum walls (diaphragm) between two neighboring tile in one ring (phi wall)
686  wallName.str("");
687  wallName << "supportWallPhi_" << iRing + 1;
688  G4Box* wall = new G4Box(wallName.str(), (wallR[iRing] - wallR[iRing - 1] - wallThick) / 2. - 1., thick / 2., wallHeight / 2.);
689  G4LogicalVolume* wallLV = new G4LogicalVolume(wall, supportMaterial, string("ARICH.") + wallName.str());
690  double r = (wallR[iRing - 1] + wallThick + wallR[iRing]) / 2.;
692 
693  // loop over layers
694  int icopyNumber = 0;
695  for (unsigned iLayer = 1; iLayer < nLayer + 1; iLayer++) {
696  //int iSlot = 1;
697  double iphi = 0;
698 
699  int iicolumn = 0;
700  // loop over phi (over tile slots from same ring)
701  while (iphi < 2 * M_PI - 0.0001) {
702 
703  // Define layer thickness as -1 in case it will not be defined
704  // below in the code with a appropriate value - Geant4 will trigger error
705  double layerThick = -1.0;
706  double tileUpThick = -1.0;
707  double tileDownThick = -1.0;
708  //cout<<" double layerThick = aeroGeo.getLayerThickness(iLayer) = "<<layerThick<<endl;
709 
710  // Define material and thickness
711  G4Material* aeroMaterial = NULL;
712  int ati_ring = iRing;
713  int ati_column = iicolumn + 1;
714  int ati_layerN = iLayer - 1;
715 
716  //cout<<setw(5)<<ati_layerN<<setw(5)<<ati_ring<<setw(5)<<ati_column<<endl;
717  if (detectorGeo.getAerogelPlane().getFullAerogelMaterialDescriptionKey() == 1) {
718  aeroMaterial = Materials::get(aeroGeo.getTileMaterialName(ati_ring, ati_column, ati_layerN).c_str());
719  } else {
720  B2ERROR("GeoARICHCreator::buildAerogelPlaneWithIndividualTilesProp --> getFullAerogelMaterialDescriptionKey() is wrong");
721  }
722  //cout<<"-----------------"<<endl
723  // <<"iLayer = "<<iLayer<<endl
724  // <<aeroMaterial->GetName()<<endl;
725  //cout<<setw(5)<<ati_layerN<<setw(5)<<ati_ring<<setw(5)<<ati_column<<endl;
726  //aeroMaterial->GetMaterialPropertiesTable()->DumpTable();
727  //zcout<<"ooooooooooooooooo"<<endl;
728  layerThick = aeroGeo.getTileThickness(ati_ring, ati_column, ati_layerN);
729  tileUpThick = aeroGeo.getTileThickness(ati_ring, ati_column, 0);
730  tileDownThick = aeroGeo.getTileThickness(ati_ring, ati_column, 1);
731 
732  // Placement of the aluminum walls (diaphragm) between two neighboring tile in one ring (phi wall)
733  // please note that we place single wall for two aerogel layers
734  G4ThreeVector trans(r * cos(iphi), r * sin(iphi), (thick - imgTubeLen) / 2.);
735  G4RotationMatrix Ra;
736  Ra.rotateZ(iphi);
737  if (iLayer == 1) {
738  new G4PVPlacement(G4Transform3D(Ra, trans), //transformation
739  wallLV, //its logical
740  string("ARICH.") + wallName.str(), //name
741  aerogelPlaneLV, //mother logical
742  false, //always false
743  icopyNumber); //should be set to 0 for the first volume of a given type.
745 
746  // In case of individual thickness
747  // (if aeroGeo.getFullAerogelMaterialDescriptionKey() == 1) or
748  // (if aeroGeo.getFullAerogelMaterialDescriptionKey() == 2)
749  // of the tiles we need to define compensation volume with ARICH air.
750  // This volume situated between aerogel tile and image plane (imgTube).
751  // This volume have same shape as earogel tile but different thickness.
752  // Build Compensation tiles only one time
753  // Compensation tile shape
754  double compTileUpThick = wallHeight - tileUpThick - tileDownThick;
755  std::stringstream compTileName;
756  //compTileName << "aerogelCompTile_" << ati_layerN << "_" << ati_ring << "_" << ati_column;
757  // In the end of the name we have Layer(L), Ring(R), Slot/Column(S) id's
758  // L : 1-2
759  // R : 1-4
760  // S : 1-22 @ R = 1
761  // S : 1-28 @ R = 2
762  // S : 1-34 @ R = 3
763  // S : 1-40 @ R = 4
764  compTileName << "aerogelCompTile_" << iLayer << "_" << ati_ring << "_" << ati_column;
765  //cout<<compTileName.str()<<endl;
766  G4Tubs* compTileShape = new G4Tubs(compTileName.str(), //name
767  wallR[iRing - 1] + wallThick + tileGap, //Rmin
768  wallR[iRing] - tileGap, //Rmax
769  compTileUpThick / 2.0, //Thikness
770  (tileGap + wallThick / 2.0) / wallR[iRing], //phi start
771  dphi - (2.0 * tileGap + wallThick) / wallR[iRing]); //delta phi
772 
773  // Logical volume of the compensation tiles
774  G4LogicalVolume* compTileLV = new G4LogicalVolume(compTileShape, //Its solid
775  imgMaterial, //G4 material
776  string("ARICH.") + compTileName.str()); //name
777 
778  // Placement of the compensation tiles
779  G4ThreeVector transCompTile(0, 0, (thick + wallHeight - compTileUpThick - imgTubeLen) / 2.0);
780  G4RotationMatrix compRa;
781  compRa.rotateZ(iphi);
782  new G4PVPlacement(G4Transform3D(compRa, transCompTile), //transformation
783  compTileLV, //its logical
784  string("ARICH.") + compTileName.str(), //name
785  aerogelPlaneLV, //mother logical
786  false, //always false
787  0); //should be set to 0 for the first volume of a given type.
788  }
789 
790  // Tile shape
791  std::stringstream tileName;
792  tileName << "aerogelTile_" << iLayer << "_" << ati_ring << "_" << ati_column;
793  //cout<<tileName.str()<<endl;
794  G4Tubs* tileShape = new G4Tubs(tileName.str(), //name
795  wallR[iRing - 1] + wallThick + tileGap, //Rmin
796  wallR[iRing] - tileGap, //Rmax
797  layerThick / 2.0, //Thikness
798  (tileGap + wallThick / 2.0) / wallR[iRing], //phi start
799  dphi - (2.0 * tileGap + wallThick) / wallR[iRing]); //delta phi
800 
801  // Logical volume of the aerogel tiles
802  G4LogicalVolume* tileLV = new G4LogicalVolume(tileShape, //Its solid
803  aeroMaterial, //G4 material
804  string("ARICH.") + tileName.str()); //name
805 
806  // Placement of the aerogel tiles
807  double zLayer = 0.0;
808  if (iLayer == 2)
809  zLayer = tileUpThick;
810  G4ThreeVector transTile(0, 0, (thick + layerThick - wallHeight - imgTubeLen) / 2.0 + zLayer);
811  new G4PVPlacement(G4Transform3D(Ra, transTile), //transformation
812  tileLV, //its logical
813  string("ARICH.") + tileName.str(), //name
814  aerogelPlaneLV, //mother logical
815  false, //always false
816  0); //should be set to 0 for the first volume of a given type.
818 
819  iphi += dphi;
820  //iSlot++;
821  icopyNumber++;
822  iicolumn++;
823  }
824  }
825  }
826 
827  // Placement of the support tube
828  new G4PVPlacement(G4Translate3D(0., 0., -(wallHeight + imgTubeLen) / 2.), //transformation
829  supportTubeLV, //its logical
830  "ARICH.AerogelSupportPlate", //name
831  aerogelPlaneLV, //mother logical
832  false, //always false
833  0); //should be set to 0 for the first volume of a given type.
834 
835  // Placement of the imaginary tube after aerogel layers (used as volume to which tracks are extrapolated by ext module)
836  new G4PVPlacement(G4Translate3D(0., 0., (wallHeight + thick) / 2.), //transformation
837  imgTubeLV, //its logical
838  "ARICH.AerogelImgPlate", //name
839  aerogelPlaneLV, //mother logical
840  false, //always false
841  0); //should be set to 0 for the first volume of a given type.
842 
843  return aerogelPlaneLV;
844 
845  }
846 
847  G4LogicalVolume* GeoARICHCreator::buildHAPD(const ARICHGeoHAPD& hapdGeo)
848  {
849 
850  // get module materials
851  string wallMat = hapdGeo.getWallMaterial();
852  string winMat = hapdGeo.getWinMaterial();
853  string apdMat = hapdGeo.getAPDMaterial();
854  string fillMat = hapdGeo.getFillMaterial();
855  string febMat = hapdGeo.getFEBMaterial();
856  G4Material* wallMaterial = Materials::get(wallMat);
857  G4Material* windowMaterial = Materials::get(winMat);
858  G4Material* apdMaterial = Materials::get(apdMat);
859  G4Material* fillMaterial = Materials::get(fillMat);
860  G4Material* febMaterial = Materials::get(febMat);
861  G4Material* moduleFill = Materials::get("ARICH_Air");
862 
863  // check that module window material has specified refractive index
864  double wref = getAvgRINDEX(windowMaterial);
865  if (!wref) B2WARNING("Material '" << winMat <<
866  "', required for ARICH photon detector window as no specified refractive index. Continuing, but no photons in ARICH will be detected.");
867 
868  // get module dimensions
869  const double hapdSizeX = hapdGeo.getSizeX();
870  const double hapdSizeY = hapdGeo.getSizeY();
871  const double hapdSizeZ = hapdGeo.getSizeZ();
872  const double wallThick = hapdGeo.getWallThickness();
873  const double winThick = hapdGeo.getWinThickness();
874  const double apdSizeX = hapdGeo.getAPDSizeX();
875  const double apdSizeY = hapdGeo.getAPDSizeY();
876  const double apdSizeZ = hapdGeo.getAPDSizeZ();
877  const double botThick = wallThick;
878  const double modHeight = hapdGeo.getModuleSizeZ();
879 
880  // module master volume
881  G4Box* moduleBox = new G4Box("moduleBox", hapdSizeX / 2., hapdSizeY / 2., modHeight / 2.);
882  G4LogicalVolume* lmoduleBox = new G4LogicalVolume(moduleBox, moduleFill, "ARICH.HAPDModule");
883 
884  // build HAPD box
885  G4Box* hapdBox = new G4Box("hapdBox", hapdSizeX / 2., hapdSizeY / 2., hapdSizeZ / 2.);
886  G4LogicalVolume* lhapdBox = new G4LogicalVolume(hapdBox, fillMaterial, "ARICH.HAPD");
887 
888  // build HAPD walls
889  G4Box* tempBox2 = new G4Box("tempBox2", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick,
890  hapdSizeZ / 2. + 0.1); // Dont't care about "+0.1", needs to be there.
891  G4SubtractionSolid* moduleWall = new G4SubtractionSolid("Box-tempBox", hapdBox, tempBox2);
892  G4LogicalVolume* lmoduleWall = new G4LogicalVolume(moduleWall, wallMaterial, "ARICH.HAPDWall");
893  setColor(*lmoduleWall, "rgb(1.0,0.0,0.0,1.0)");
894  new G4PVPlacement(G4Transform3D(), lmoduleWall, "ARICH.HAPDWall", lhapdBox, false, 1);
895 
896  // build HAPD window
897  G4Box* winBox = new G4Box("winBox", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick, winThick / 2.);
898  G4LogicalVolume* lmoduleWin = new G4LogicalVolume(winBox, windowMaterial, "ARICH.HAPDWindow");
899  setColor(*lmoduleWin, "rgb(0.7,0.7,0.7,1.0)");
900  lmoduleWin->SetSensitiveDetector(m_sensitive);
901  G4Transform3D transform = G4Translate3D(0., 0., (-hapdSizeZ + winThick) / 2.);
902  new G4PVPlacement(transform, lmoduleWin, "ARICH.HAPDWindow", lhapdBox, false, 1);
903 
904  // build module bottom
905  G4Box* botBox = new G4Box("botBox", hapdSizeX / 2. - wallThick, hapdSizeY / 2. - wallThick, botThick / 2.);
906  G4LogicalVolume* lmoduleBot = new G4LogicalVolume(botBox, wallMaterial, "ARICH.HAPDBottom");
907  setColor(*lmoduleBot, "rgb(0.0,1.0,0.0,1.0)");
908  G4Transform3D transform1 = G4Translate3D(0., 0., (hapdSizeZ - botThick) / 2.);
909  new G4PVPlacement(transform1, lmoduleBot, "ARICH.HAPDBottom", lhapdBox, false, 1);
910 
911  // build apd
912  G4Box* apdBox = new G4Box("apdBox", apdSizeX / 2., apdSizeY / 2., apdSizeZ / 2.);
913  G4LogicalVolume* lApd = new G4LogicalVolume(apdBox, apdMaterial, "ARICH.HAPDApd");
914  if (m_isBeamBkgStudy) lApd->SetSensitiveDetector(new BkgSensitiveDetector("ARICH", 1));
915 
916  // add APD surface optical properties
917  Materials& materials = Materials::getInstance();
918 
919  G4OpticalSurface* optSurf = materials.createOpticalSurface(hapdGeo.getAPDSurface());
920 
921  new G4LogicalSkinSurface("apdSurface", lApd, optSurf);
922  G4Transform3D transform2 = G4Translate3D(0., 0., (hapdSizeZ - apdSizeZ) / 2. - botThick);
923  new G4PVPlacement(transform2, lApd, "ARICH.HAPDApd", lhapdBox, false, 1);
924 
925  // build FEB
926  double febSizeX = hapdGeo.getFEBSizeX();
927  double febSizeY = hapdGeo.getFEBSizeY();
928  double febSizeZ = hapdGeo.getFEBSizeZ();
929  G4Box* febBox = new G4Box("febBox", febSizeX / 2., febSizeY / 2., febSizeZ / 2.);
930  G4LogicalVolume* lfeb = new G4LogicalVolume(febBox, febMaterial, "ARICH.HAPDFeb");
931  if (m_isBeamBkgStudy) lfeb->SetSensitiveDetector(new BkgSensitiveDetector("ARICH"));
932  setColor(*lfeb, "rgb(0.0,0.6,0.0,1.0)");
933  G4Transform3D transform3 = G4Translate3D(0., 0., (modHeight - febSizeZ) / 2.);
934  new G4PVPlacement(transform3, lfeb, "ARICH.HAPDFeb", lmoduleBox, false, 1);
935  G4Transform3D transform4 = G4Translate3D(0., 0., - (modHeight - hapdSizeZ) / 2.);
936  new G4PVPlacement(transform4, lhapdBox, "ARICH.HAPD", lmoduleBox, false, 1);
937 
938  return lmoduleBox;
939 
940  }
941 
942  G4LogicalVolume* GeoARICHCreator::buildDetectorPlane(const ARICHGeometryConfig& detectorGeo)
943  {
944 
945  const ARICHGeoHAPD& hapdGeo = detectorGeo.getHAPDGeometry();
946  G4LogicalVolume* hapdLV = buildHAPD(hapdGeo);
947 
948  const ARICHGeoDetectorPlane& detGeo = detectorGeo.getDetectorPlane();
949 
950  G4Tubs* detTube = new G4Tubs("detTube", detGeo.getRingR(1) - hapdGeo.getSizeX() * 1.4 / 2.,
951  detGeo.getRingR(detGeo.getNRings()) + hapdGeo.getSizeX() * 1.4 / 2. , hapdGeo.getModuleSizeZ() / 2., 0, 2 * M_PI);
952  G4LogicalVolume* detPlaneLV = new G4LogicalVolume(detTube, Materials::get("ARICH_Air"), "ARICH.detectorPlane");
953 
954  unsigned nSlots = detGeo.getNSlots();
955 
956  for (unsigned iSlot = 1; iSlot < nSlots + 1; iSlot++) {
957  if (!m_modInfo->isInstalled(iSlot)) continue;
958  double r = detGeo.getSlotR(iSlot);
959  double phi = detGeo.getSlotPhi(iSlot);
960  G4ThreeVector trans(r * cos(phi), r * sin(phi), 0);
961  G4RotationMatrix Ra;
962  Ra.rotateZ(phi);
963  G4ThreeVector trans1(r * cos(phi), r * sin(phi), 0.0);
964  new G4PVPlacement(G4Transform3D(Ra, trans1), hapdLV, "ARICH.HAPDModule", detPlaneLV, false, iSlot);
965  }
966 
967  return detPlaneLV;
968 
969  }
970 
971  G4LogicalVolume* GeoARICHCreator::buildMerger(const ARICHGeoMerger& mergerGeo)
972  {
973 
974  // This is a screw hole on the merger board.
975  // This high precision of the geometry description is needed
976  // only for the correct placement of the merger cooling bodies.
977  // Volume to subtract (screw hole on the merger board)
978  double screwholeR = mergerGeo.getMergerPCBscrewholeR() * mm;
979  double screwholedY = mergerGeo.getMergerPCBscrewholePosdY() * mm;
980  double screwholedX1 = mergerGeo.getMergerPCBscrewholePosdX1() * mm;
981  double screwholedX2 = mergerGeo.getMergerPCBscrewholePosdX2() * mm;
982 
983  G4VSolid* screwHoleTubeSubtracted_solid = new G4Tubs("screwHoleTubeSubtracted_solid",
984  0.0,
985  screwholeR,
986  mergerGeo.getThickness() * mm / 2.0,
987  0, 360.0 * deg);
988 
989  // Volume to add (merger box)
990  G4Box* merger_solid = new G4Box("merger_solid",
991  mergerGeo.getSizeW() * mm / 2.0,
992  mergerGeo.getSizeL() * mm / 2.0,
993  mergerGeo.getThickness() * mm / 2.0);
994 
995  G4RotationMatrix Ra_sub;
996  G4ThreeVector Ta_sub;
997  G4Transform3D Tr_sub;
998  Ta_sub.setX(-mergerGeo.getSizeW() * mm / 2.0 + screwholedX1);
999  Ta_sub.setY(mergerGeo.getSizeL() * mm / 2.0 - screwholedY);
1000  Ta_sub.setZ(0.0);
1001  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1002  G4SubtractionSolid* substraction_solid = new G4SubtractionSolid("substraction_solid", merger_solid, screwHoleTubeSubtracted_solid,
1003  Tr_sub);
1004  Ta_sub.setX(mergerGeo.getSizeW() * mm / 2.0 - screwholedX2);
1005  Ta_sub.setY(mergerGeo.getSizeL() * mm / 2.0 - screwholedY);
1006  Ta_sub.setZ(0.0);
1007  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1008  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, screwHoleTubeSubtracted_solid, 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 + screwholedX1);
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 
1020  return new G4LogicalVolume(substraction_solid, Materials::get(mergerGeo.getMergerPCBMaterialName()), "ARICH.mergerPCB");
1021  }
1022 
1023 
1024  G4LogicalVolume* GeoARICHCreator::buildMergerCooling(unsigned iType)
1025  {
1026 
1027  if (!m_mergerCooling) {
1028  B2WARNING("ARICH geometry: no data available for merger " << iType << " cooling body geometry. Cooling body will not be placed.");
1029  return NULL;
1030  }
1031 
1032  std::stringstream shpName;
1033  shpName << "TessellatedSolid_" << + iType;
1034 
1035  G4TessellatedSolid* volume_solid = new G4TessellatedSolid(shpName.str().c_str());
1036 
1037  G4ThreeVector point_1;
1038  G4ThreeVector point_2;
1039  G4ThreeVector point_3;
1040 
1041  tessellatedSolidStr mergerCoolingStr = m_mergerCooling->getMergerCoolingBodiesInfo(iType);
1042 
1043  if (mergerCoolingStr.nCells == 0) {
1044  B2WARNING("ARICH geometry: no data available for merger " << iType << " cooling body geometry. Cooling body will not be placed.");
1045  return NULL;
1046  }
1047 
1048  for (unsigned int i = 0; i < mergerCoolingStr.nCells; i++) {
1049  //
1050  point_1.setX(mergerCoolingStr.posV1[0][i]);
1051  point_1.setY(mergerCoolingStr.posV1[1][i]);
1052  point_1.setZ(mergerCoolingStr.posV1[2][i]);
1053  //
1054  point_2.setX(mergerCoolingStr.posV2[0][i]);
1055  point_2.setY(mergerCoolingStr.posV2[1][i]);
1056  point_2.setZ(mergerCoolingStr.posV2[2][i]);
1057  //
1058  point_3.setX(mergerCoolingStr.posV3[0][i]);
1059  point_3.setY(mergerCoolingStr.posV3[1][i]);
1060  point_3.setZ(mergerCoolingStr.posV3[2][i]);
1061  //
1062  G4TriangularFacet* facet = new G4TriangularFacet(point_1, point_2, point_3, ABSOLUTE);
1063  volume_solid->AddFacet((G4VFacet*) facet);
1064  }
1065 
1066  volume_solid->SetSolidClosed(true);
1067  std::stringstream volName;
1068  volName << "ARICH.mergerCooling_" << + iType;
1069  G4LogicalVolume* volume_logical = new G4LogicalVolume(volume_solid,
1070  Materials::get(m_mergerCooling->getMergerCoolingBodiesMaterialName()), volName.str().c_str());
1071 
1072  setColor(*volume_logical, "rgb(0.6,0.0,0.2,1.0)"); //From the top (farther from IP, downstream)
1073 
1074  return volume_logical;
1075  }
1076 
1077  G4LogicalVolume* GeoARICHCreator::buildMergerEnvelope(const ARICHGeoMerger& mergerGeo, int iType)
1078  {
1079  // Volume to add single merger and merger cooling body envelope box
1080  G4Box* singlemergerenvelope_solid = new G4Box("singlemergerenvelope_solid",
1081  mergerGeo.getSingleMergerEnvelopeSizeW() * mm / 2.0,
1082  mergerGeo.getSingleMergerEnvelopeSizeL() * mm / 2.0,
1083  mergerGeo.getSingleMergerEnvelopeThickness() * mm / 2.0);
1084  std::stringstream volName;
1085  volName << "ARICH.singleMergerEnvelope_" << + iType;
1086  return new G4LogicalVolume(singlemergerenvelope_solid, Materials::get("ARICH_Air"), volName.str().c_str());
1087  }
1088 
1090  {
1091 
1092  const ARICHGeoMerger& mergerGeo = detectorGeo.getMergerGeometry();
1093 
1094  if (mergerGeo.getSingleMergerEnvelopeSizeW() < 1e-9) {
1095  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).");
1096  return NULL;
1097  }
1098 
1099  G4Tubs* envelope_solid = new G4Tubs("envelope_solid", mergerGeo.getEnvelopeInnerRadius() * mm,
1100  mergerGeo.getEnvelopeOuterRadius() * mm, mergerGeo.getEnvelopeThickness() * mm / 2.0, 0.0, 2.0 * M_PI);
1101 
1102  G4LogicalVolume* envelope_logical = new G4LogicalVolume(envelope_solid, Materials::get("ARICH_Air"), "ARICH.mergerEnvelope");
1103 
1104 
1105  G4LogicalVolume* merger_logical = buildMerger(mergerGeo);
1106  G4LogicalVolume* mergerCooling_logical[12] = {NULL};
1107  G4LogicalVolume* mergerEnvelope_logical[12] = {NULL};
1108  G4ThreeVector TaPCB(mergerGeo.getSingleMergeEnvelopePosition().X() * mm, mergerGeo.getSingleMergeEnvelopePosition().Y() * mm,
1109  mergerGeo.getSingleMergeEnvelopePosition().Z() * mm);
1110  G4RotationMatrix RaPCB;
1111  G4ThreeVector TaMergerCooling(mergerGeo.getSingleMergeEnvelopePosition().X() * mm,
1112  mergerGeo.getSingleMergeEnvelopePosition().Y() * mm, -mergerGeo.getSingleMergeEnvelopePosition().Z() * mm);
1113  G4RotationMatrix RaMergerCooling;
1114  RaMergerCooling.rotateY(180 * deg);
1115  RaMergerCooling.rotateZ(-90 * deg);
1116 
1117  // build 12 different merger+cooling body volumes
1118  for (int iType = 1; iType < 13; iType++) {
1119  mergerCooling_logical[iType - 1] = buildMergerCooling(iType);
1120  mergerEnvelope_logical[iType - 1 ] = buildMergerEnvelope(mergerGeo, iType); //Single merger and merger cooling body envelope box
1121  setColor(*mergerEnvelope_logical[iType - 1], "rgb(0.0,0.0,1.0,1.0)");
1122 
1123  new G4PVPlacement(G4Transform3D(RaPCB, TaPCB), //Transformation
1124  merger_logical, //its logical volume
1125  "ARICH.mergerPCB", //its name
1126  mergerEnvelope_logical[iType - 1], //its mother volume
1127  false, //no boolean operation
1128  iType); //copy number
1129 
1130  if (mergerCooling_logical[iType - 1] == NULL) continue;
1131 
1132  new G4PVPlacement(G4Transform3D(RaMergerCooling, TaMergerCooling), //Transformation
1133  mergerCooling_logical[iType - 1], //its logical volume
1134  "ARICH.mergerCooling", //its name
1135  mergerEnvelope_logical[iType - 1], //its mother volume
1136  false, //no boolean operation
1137  iType); //copy number
1138  }
1139 
1140  // place all 72 merger+cooling body packages
1141  for (unsigned iSlot = 0; iSlot < mergerGeo.getMergerSlotID().size(); iSlot++) {
1142 
1143  int type = 1; // if no merger cooling is available...
1144  if (m_mergerCooling) type = (int)m_mergerCooling->getMergerCoolingPositionID().at(iSlot);
1145 
1146  //Placement of the single volume envelope
1147  G4ThreeVector Ta(mergerGeo.getMergerPosR().at(iSlot) * mm * cos(mergerGeo.getMergerAngle().at(iSlot) * deg),
1148  mergerGeo.getMergerPosR().at(iSlot) * mm * sin(mergerGeo.getMergerAngle().at(iSlot) * deg),
1149  mergerGeo.getSingleMergerenvelopeDeltaZ().at(iSlot) * mm);
1150  G4RotationMatrix Ra;
1151  Ra.rotateZ(mergerGeo.getMergerAngle().at(iSlot) * deg + mergerGeo.getMergerOrientation().at(iSlot) * deg);
1152  new G4PVPlacement(G4Transform3D(Ra, Ta), //Transformation
1153  mergerEnvelope_logical[type - 1], //its logical volume
1154  "ARICH.singleMergerEnvelope", //its name
1155  envelope_logical, //its mother volume
1156  false, //no boolean operation
1157  iSlot); //copy number
1158  }
1159 
1160  return envelope_logical;
1161 
1162  }
1163 
1164  G4LogicalVolume* GeoARICHCreator::buildCables(const ARICHGeoCablesEnvelope& cablesGeo)
1165  {
1166 
1167  G4Tubs* cablesEnvelope_solid = new G4Tubs("cablesEnvelope_solid",
1168  cablesGeo.getEnvelopeInnerRadius(),
1169  cablesGeo.getEnvelopeOuterRadius(),
1170  cablesGeo.getEnvelopeThickness() / 2.0,
1171  0.0,
1172  2.0 * M_PI);
1173  G4LogicalVolume* cablesEnvelope_logical = new G4LogicalVolume(cablesEnvelope_solid,
1175  "ARICH.cablesEnvelope");
1176 
1177  return cablesEnvelope_logical;
1178 
1179  }
1180 
1181  G4LogicalVolume* GeoARICHCreator::buildFEBCoolingBody(const ARICHGeoFEBCooling& coolingv2Geo)
1182  {
1183 
1184  //FEB aluminum cooling envelope for single object
1185  double feb_alcooling_singleObjectEnvelope_sizeX = (2 * coolingv2Geo.getSmallSquareSize() + coolingv2Geo.getBigSquareSize() + 2.0) *
1186  mm;
1187  double feb_alcooling_singleObjectEnvelope_sizeY = feb_alcooling_singleObjectEnvelope_sizeX * mm;
1188  double feb_alcooling_singleObjectEnvelope_sizeZ = (coolingv2Geo.getSmallSquareThickness() + coolingv2Geo.getRectangleThickness()) *
1189  mm;
1190 
1191  double feb_alcooling_box1_sizeX = coolingv2Geo.getSmallSquareSize() * mm;
1192  double feb_alcooling_box1_sizeY = feb_alcooling_box1_sizeX;
1193  double feb_alcooling_box1_sizeZ = coolingv2Geo.getSmallSquareThickness() * mm;
1194 
1195  double feb_alcooling_box2_sizeX = coolingv2Geo.getBigSquareSize() * mm;
1196  double feb_alcooling_box2_sizeY = feb_alcooling_box2_sizeX;
1197  double feb_alcooling_box2_sizeZ = coolingv2Geo.getBigSquareThickness() * mm;
1198 
1199  double feb_alcooling_box3_sizeX = coolingv2Geo.getRectangleW() * mm;
1200  double feb_alcooling_box3_sizeY = coolingv2Geo.getRectangleL() * mm;
1201  double feb_alcooling_box3_sizeZ = coolingv2Geo.getRectangleThickness() * mm;
1202 
1203  double feb_alcooling_box1_X0 = feb_alcooling_box2_sizeX / 2.0 + feb_alcooling_box1_sizeX / 2.0;
1204  double feb_alcooling_box1_Y0 = feb_alcooling_box2_sizeY / 2.0 + feb_alcooling_box1_sizeY / 2.0;
1205  double feb_alcooling_box1_Z0 = 0.0 * mm;
1206 
1207  //double feb_alcooling_box2_X0 = 0.0*mm;
1208  //double feb_alcooling_box2_Y0 = 0.0*mm;
1209  //double feb_alcooling_box2_Z0 = 0.0*mm;
1210 
1211  double feb_alcooling_box3_X0 = coolingv2Geo.getRectangleDistanceFromCenter() / sqrt(2.0) * mm;
1212  double feb_alcooling_box3_Y0 = feb_alcooling_box3_X0;
1213  double feb_alcooling_box3_Z0 = feb_alcooling_box1_sizeZ / 2.0 + feb_alcooling_box3_sizeZ / 2.0;
1214  double feb_alcooling_box3_angle = 45.0 * deg;
1215 
1216  G4RotationMatrix Ra;
1217  G4ThreeVector Ta;
1218  G4Transform3D Tr;
1219 
1220  //
1221  // Define single FEB aluminum cooling envelope
1222  //
1223  G4VSolid* feb_alcoolingEnvelope_solid = new G4Box("feb_alcoolingEnvelope_solid",
1224  feb_alcooling_singleObjectEnvelope_sizeX / 2.0,
1225  feb_alcooling_singleObjectEnvelope_sizeY / 2.0,
1226  feb_alcooling_singleObjectEnvelope_sizeZ / 2.0);
1227  G4LogicalVolume* feb_alcoolingEnvelope_logical = new G4LogicalVolume(feb_alcoolingEnvelope_solid, Materials::get("Air"),
1228  "feb_alcoolingEnvelope_logical");
1229 
1230  G4VSolid* feb_alcooling_box1_solid = new G4Box("feb_alcooling_box1_solid", feb_alcooling_box1_sizeX / 2.0,
1231  feb_alcooling_box1_sizeY / 2.0, feb_alcooling_box1_sizeZ / 2.0);
1232  G4VSolid* feb_alcooling_box2_solid = new G4Box("feb_alcooling_box2_solid", feb_alcooling_box2_sizeX / 2.0,
1233  feb_alcooling_box2_sizeY / 2.0, feb_alcooling_box2_sizeZ / 2.0);
1234  G4VSolid* feb_alcooling_box3_solid = new G4Box("feb_alcooling_box3_solid", feb_alcooling_box3_sizeX / 2.0,
1235  feb_alcooling_box3_sizeY / 2.0, feb_alcooling_box3_sizeZ / 2.0);
1236 
1237  //
1238  //Box 1 A
1239  //
1240  Ta.setX(feb_alcooling_box1_X0);
1241  Ta.setY(feb_alcooling_box1_Y0);
1242  Ta.setZ(feb_alcooling_box1_Z0);
1243  Tr = G4Transform3D(Ra, Ta);
1244  G4UnionSolid* feb_alcooling_assembly01_solid = new G4UnionSolid("feb_alcooling_assembly01_solid", feb_alcooling_box2_solid,
1245  feb_alcooling_box1_solid, Tr);
1246  //
1247  //Box 1 B
1248  //
1249  Ta.setX(-feb_alcooling_box1_X0);
1250  Ta.setY(-feb_alcooling_box1_Y0);
1251  Ta.setZ(feb_alcooling_box1_Z0);
1252  Tr = G4Transform3D(Ra, Ta);
1253  G4UnionSolid* feb_alcooling_assembly02_solid = new G4UnionSolid("feb_alcooling_assembly02_solid", feb_alcooling_assembly01_solid,
1254  feb_alcooling_box1_solid, Tr);
1255  //
1256  //Box 3 A
1257  //
1258  Ta.setX(feb_alcooling_box3_X0);
1259  Ta.setY(feb_alcooling_box3_Y0);
1260  Ta.setZ(feb_alcooling_box3_Z0);
1261  Ra.rotateZ(-feb_alcooling_box3_angle);
1262  Tr = G4Transform3D(Ra, Ta);
1263  G4UnionSolid* feb_alcooling_assembly03_solid = new G4UnionSolid("feb_alcooling_assembly03_solid", feb_alcooling_assembly02_solid,
1264  feb_alcooling_box3_solid, Tr);
1265  Ra.rotateZ(feb_alcooling_box3_angle);
1266  //
1267  //Box 3 B
1268  //
1269  Ta.setX(-feb_alcooling_box3_X0);
1270  Ta.setY(-feb_alcooling_box3_Y0);
1271  Ta.setZ(feb_alcooling_box3_Z0);
1272  Ra.rotateZ(-feb_alcooling_box3_angle);
1273  Tr = G4Transform3D(Ra, Ta);
1274  G4UnionSolid* feb_alcooling_assembly_solid = new G4UnionSolid("feb_alcooling_assembly_solid", feb_alcooling_assembly03_solid,
1275  feb_alcooling_box3_solid, Tr);
1276  Ra.rotateZ(feb_alcooling_box3_angle);
1277 
1278  G4LogicalVolume* feb_alcooling_assembly_logical = new G4LogicalVolume(feb_alcooling_assembly_solid, Materials::get("Al"),
1279  "feb_alcooling_assembly_logical");
1280  Ta.setX(0.0);
1281  Ta.setY(0.0);
1282  Ta.setZ(-feb_alcooling_box3_sizeZ / 2.0);
1283 
1284  Tr = G4Transform3D(Ra, Ta);
1285  new G4PVPlacement(Tr, //Transformation
1286  feb_alcooling_assembly_logical, //its logical volume
1287  "feb_alcooling_assembly", //its name
1288  feb_alcoolingEnvelope_logical, //its mother volume
1289  false, //no boolean operation
1290  0); //copy number
1291 
1292  return feb_alcoolingEnvelope_logical;
1293  }
1294 
1295  G4LogicalVolume* GeoARICHCreator::buildCoolingTube(const unsigned i_volumeID, const ARICHGeoCooling& coolingGeo)
1296  {
1297 
1298  B2ASSERT("ARICH cooling geometry ID (G4Tube) is wrong : coolingGeo.getCoolingGeometryID.at(i_volumeID) != 1",
1299  coolingGeo.getCoolingGeometryID().at(i_volumeID) == 1);
1300  G4Tubs* coolingTube_solid = new G4Tubs("coolingTube_solid",
1301  coolingGeo.getRmin() * mm,
1302  coolingGeo.getRmax() * mm,
1303  coolingGeo.getCoolingL().at(i_volumeID) * mm / 2.0,
1304  0.0,
1305  2.0 * M_PI);
1306  return new G4LogicalVolume(coolingTube_solid, Materials::get(coolingGeo.getCoolingPipeMaterialName()), "ARICH.coolingTube");
1307 
1308  }
1309 
1310  G4LogicalVolume* GeoARICHCreator::buildCoolingTorus(const unsigned i_volumeID, const ARICHGeoCooling& coolingGeo)
1311  {
1312 
1313  B2ASSERT("ARICH cooling geometry ID (G4Torus) is wrong : coolingGeo.getCoolingGeometryID.at(i_volumeID) != 2",
1314  coolingGeo.getCoolingGeometryID().at(i_volumeID) == 2);
1315 
1316  double pSPhi = coolingGeo.getCoolingPosPhi().at(i_volumeID) * deg - coolingGeo.getCoolingL().at(
1317  i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID) / 2.0;
1318  //double pDPhi = coolingGeo.getCoolingPosPhi().at(i_volumeID)*deg + coolingGeo.getCoolingL().at(i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID) / 2.0;
1319  double pDPhi = coolingGeo.getCoolingL().at(i_volumeID) / coolingGeo.getCoolingPosR().at(i_volumeID);
1320 
1321  G4Torus* coolingTorus_solid = new G4Torus("coolingTorus_solid", // name
1322  coolingGeo.getRmin(), // pRmin
1323  coolingGeo.getRmax(), // pRmax
1324  coolingGeo.getCoolingPosR().at(i_volumeID), // pRtor
1325  pSPhi, // pSPhi
1326  pDPhi); // pDPhi
1327 
1328  return new G4LogicalVolume(coolingTorus_solid, Materials::get(coolingGeo.getCoolingPipeMaterialName()), "ARICH.coolingTorus");
1329 
1330  }
1331 
1333  {
1334 
1335  G4Tubs* coolingEnvelope_solid = new G4Tubs("coolingEnvelope_solid",
1336  coolingGeo.getEnvelopeInnerRadius(),
1337  coolingGeo.getEnvelopeOuterRadius(),
1338  coolingGeo.getEnvelopeThickness() / 2.0,
1339  0.0,
1340  2.0 * M_PI);
1341  G4LogicalVolume* coolingEnvelope_logical = new G4LogicalVolume(coolingEnvelope_solid,
1342  Materials::get("Air"),
1343  "ARICH.coolingEnvelope");
1344 
1345  unsigned nComponents = coolingGeo.getCoolingGeometryID().size();
1346 
1347  for (unsigned i = 0; i < nComponents; i++) {
1348  double r = coolingGeo.getCoolingPosR().at(i) * mm;
1349  double phi = coolingGeo.getCoolingPosPhi().at(i) * deg;
1350  G4ThreeVector Ta;
1351  G4RotationMatrix Ra;
1352  G4LogicalVolume* coolingComponentLV;
1353  if (coolingGeo.getCoolingGeometryID().at(i) == 1) {
1354  //<!-- 1 -> G4Tubs --> build a tube
1355  Ta.set(r * cos(phi), r * sin(phi), 0);
1356  //First need to rotate around y axis
1357  //to make the tube parallel with x axis
1358  Ra.rotateY(90.0 * deg);
1359  Ra.rotateZ(coolingGeo.getCoolinRotationAngle().at(i) * deg);
1360  coolingComponentLV = buildCoolingTube(i, coolingGeo);
1361  } else if (coolingGeo.getCoolingGeometryID().at(i) == 2) {
1362  //<!-- 2 -> G4Torus --> build a torus
1363  coolingComponentLV = buildCoolingTorus(i, coolingGeo);
1364  } else {
1365  B2FATAL("ARICH cooling geometry component ID is wrong");
1366  }
1367  new G4PVPlacement(G4Transform3D(Ra, Ta), //Transformation
1368  coolingComponentLV, //its logical volume
1369  "ARICH.cooling", //its name
1370  coolingEnvelope_logical, //its mother volume
1371  false, //no boolean operation
1372  i); //copy number
1373  }
1374 
1375  return coolingEnvelope_logical;
1376 
1377  }
1378 
1379  G4LogicalVolume* GeoARICHCreator::buildCoolingTestPlate(const ARICHGeoCooling& coolingGeo)
1380  {
1381 
1382  G4Box* coolingTestPlateEnvelop_solid = new G4Box("coolingTestPlateEnvelop_solid",
1383  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1384  coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 * mm,
1385  coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 * mm);
1386  G4LogicalVolume* coolingTestPlateEnvelop_logical = new G4LogicalVolume(coolingTestPlateEnvelop_solid, Materials::get("Air"),
1387  "ARICH.coolingTestPlateEnvelop");
1388 
1389  G4Box* coolingTestPlate_solid = new G4Box("coolingTestPlate_solid",
1390  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1391  coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 * mm,
1392  coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 * mm);
1393 
1394  // Volume to subtract
1395  G4VSolid* coldTubeSubtracted_solid = new G4Tubs("coldTubeSubtracted_solid",
1396  0.0,
1397  coolingGeo.getColdTubeSubtractedR() * mm,
1398  (coolingGeo.getCoolingTestPlateslengths().X() + 1) / 2.0 * mm,
1399  0, 360.0 * deg);
1400 
1401  // Volume to add (cold tube)
1402  G4VSolid* coldTube_solid = new G4Tubs("coldTube_solid",
1403  (coolingGeo.getColdTubeR() - coolingGeo.getColdTubeWallThickness()) * mm,
1404  coolingGeo.getColdTubeR() * mm,
1405  coolingGeo.getCoolingTestPlateslengths().X() / 2.0 * mm,
1406  0, 360.0 * deg);
1407  G4LogicalVolume* coldTube_logical = new G4LogicalVolume(coldTube_solid, Materials::get(coolingGeo.getColdTubeMaterialName()),
1408  "ARICH.coldTube");
1409 
1410  G4RotationMatrix Ra_sub;
1411  G4ThreeVector Ta_sub;
1412  G4Transform3D Tr_sub;
1413  Ta_sub.setX(0.0);
1414  Ta_sub.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing()) * mm);
1415  Ta_sub.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1416  Ra_sub.rotateY(90.0 * deg);
1417  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1418  G4SubtractionSolid* substraction_solid = new G4SubtractionSolid("substraction_solid", coolingTestPlate_solid,
1419  coldTubeSubtracted_solid, Tr_sub);
1420  for (int i = 1; i < coolingGeo.getColdTubeNumber(); i++) {
1421  Ta_sub.setX(0.0);
1422  Ta_sub.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing() - coolingGeo.getColdTubeSpacing()
1423  * 2 * i) * mm);
1424  Ta_sub.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1425  Tr_sub = G4Transform3D(Ra_sub, Ta_sub);
1426  substraction_solid = new G4SubtractionSolid("substraction_solid", substraction_solid, coldTubeSubtracted_solid, Tr_sub);
1427  }
1428 
1429  G4LogicalVolume* coolingTestPlate_logical = new G4LogicalVolume(substraction_solid,
1430  Materials::get(coolingGeo.getCoolingTestPlateMaterialName()), "ARICH.coolingTestPlate");
1431 
1432  new G4PVPlacement(G4Transform3D(), //Transformation
1433  coolingTestPlate_logical, //its logical volume
1434  "ARICH.coolingTestPlate", //its name
1435  coolingTestPlateEnvelop_logical, //its mother volume
1436  false, //no boolean operation
1437  0); //copy number
1438  //Add cold tubes
1439  G4RotationMatrix Ra;
1440  G4ThreeVector Ta;
1441  G4Transform3D Tr;
1442  Ra.rotateY(90.0 * deg);
1443  for (int i = 0; i < coolingGeo.getColdTubeNumber(); i++) {
1444  Ta.setX(0.0);
1445  Ta.setY((coolingGeo.getCoolingTestPlateslengths().Y() / 2.0 - coolingGeo.getColdTubeSpacing() - coolingGeo.getColdTubeSpacing() *
1446  2 * i) * mm);
1447  Ta.setZ((coolingGeo.getCoolingTestPlateslengths().Z() / 2.0 - coolingGeo.getDepthColdTubeInPlate()) * mm);
1448  Tr = G4Transform3D(Ra, Ta);
1449  new G4PVPlacement(Tr, //Transformation
1450  coldTube_logical, //its logical volume
1451  "ARICH.coldTube", //its name
1452  coolingTestPlateEnvelop_logical, //its mother volume
1453  false, //no boolean operation
1454  0); //copy number
1455  }
1456 
1457  return coolingTestPlateEnvelop_logical;
1458 
1459  }
1460 
1462  {
1463 
1464  const ARICHGeoDetectorPlane& detGeo = detectorGeo.getDetectorPlane();
1465 
1466  G4Tubs* supportTube = new G4Tubs("supportTube", detGeo.getSupportInnerR(), detGeo.getSupportOuterR(),
1467  (detGeo.getSupportThickness() + detGeo.getSupportBackWallHeight()) / 2., 0, 2 * M_PI);
1468  G4Material* supportMaterial = Materials::get(detGeo.getSupportMaterial());
1469  G4LogicalVolume* detSupportLV = new G4LogicalVolume(supportTube, supportMaterial, "ARICH.detectorSupportPlate");
1470 
1471  G4Tubs* supportPlate = new G4Tubs("supportPlate", detGeo.getSupportInnerR(), detGeo.getSupportOuterR(),
1472  detGeo.getSupportThickness() / 2., 0, 2 * M_PI);
1473 
1474  G4Box* hole = new G4Box("hole", detGeo.getModuleHoleSize() / 2., detGeo.getModuleHoleSize() / 2.,
1475  detGeo.getSupportThickness() / 2.); // +1 for thickness for subtraction solid
1476  G4LogicalVolume* holeLV = new G4LogicalVolume(hole, Materials::get("Air"), "ARICH.detectorSupportHole");
1477 
1478  int nRings = detGeo.getNRings();
1479  std::vector<G4LogicalVolume*> hapdBackRadialWallLV;
1480  double backWallThick = detGeo.getSupportBackWallThickness();
1481  double backWallHeight = detGeo.getSupportBackWallHeight();
1482 
1483  G4LogicalVolume* hapdSupportPlateLV = new G4LogicalVolume(supportPlate, supportMaterial, "hapdSupport");
1484 
1485  std::vector<double> wallR;
1486  wallR.assign(nRings + 1, 0);
1487  std::vector<double> thickR;
1488  thickR.assign(nRings + 1, 0);
1489 
1490  for (int i = 1; i < nRings; i++) {
1491  double rm1 = detGeo.getRingR(i);
1492  double rp1 = detGeo.getRingR(i + 1);
1493  wallR[i] = (rp1 + rm1) / 2.;
1494  if (i == 1) {
1495  wallR[0] = rm1 - (rp1 - rm1) / 2.;
1496  }
1497  if (i == nRings - 1) {
1498  wallR[i + 1] = rp1 + (rp1 - rm1) / 2.;
1499  }
1500  }
1501 
1502  for (int i = 0; i < nRings + 1; i++) {
1503  std::stringstream ringName1;
1504  ringName1 << "backWall_" << i;
1505  thickR[i] = backWallThick;
1506  if (i == nRings) {
1507  thickR[i] = 2.*backWallThick;
1508  wallR[i] = detGeo.getSupportOuterR() - thickR[i] / 2.;
1509  }
1510  G4Tubs* backTube = new G4Tubs("hapdBackRing", wallR[i] - thickR[i] / 2., wallR[i] + thickR[i] / 2., backWallHeight / 2., 0,
1511  2 * M_PI);
1512  G4LogicalVolume* hapdBackTubeLV = new G4LogicalVolume(backTube, supportMaterial, "backTube");
1513  G4Transform3D transform3 = G4Translate3D(0., 0., detGeo.getSupportThickness() / 2.);
1514  new G4PVPlacement(transform3, hapdBackTubeLV, "backTube", detSupportLV, false, 1);
1515  if (i == 0) continue;
1516 
1517  G4Box* backRadial = new G4Box("backRadialBox", (wallR[i] - wallR[i - 1] - thickR[i] / 2. - thickR[i - 1] / 2.) / 2. - 1.,
1518  backWallThick / 2., backWallHeight / 2.);
1519  hapdBackRadialWallLV.push_back(new G4LogicalVolume(backRadial, supportMaterial, ringName1.str().c_str()));
1520  }
1521 
1522  G4SubtractionSolid* substraction = NULL;
1523  unsigned nSlots = detGeo.getNSlots();
1524 
1525  if (nSlots > detectorGeo.getFEBCoolingGeometry().getFebcoolingv2GeometryID().size()) {
1526  B2WARNING("GeoARICHCreator: No FEB colling body geometry available so they will not be placed (ARICHGeometryConfig with ClasDef>4 is needed).");
1527  return detSupportLV;
1528  }
1529 
1530  // drill holes in support plate
1531  // add FEB cooling bodies
1532  for (unsigned iSlot = 1; iSlot < nSlots + 1; iSlot++) {
1533  unsigned iRing = detGeo.getSlotRing(iSlot);
1534  double r = (wallR[iRing] + wallR[iRing - 1]) / 2. - (thickR[iRing] - thickR[iRing - 1]) / 2.;
1535 
1536  double phi = detGeo.getSlotPhi(iSlot);
1537  G4ThreeVector trans(r * cos(phi), r * sin(phi), 0);
1538  G4RotationMatrix Ra;
1539  Ra.rotateZ(phi);
1540 
1541  new G4PVPlacement(G4Transform3D(Ra, trans), holeLV, "hole", hapdSupportPlateLV, false, iSlot);
1542  if (substraction) substraction = new G4SubtractionSolid("Box+CylinderMoved", substraction, hole, G4Transform3D(Ra, trans));
1543  else substraction = new G4SubtractionSolid("Box+CylinderMoved", supportPlate, hole, G4Transform3D(Ra, trans));
1544 
1545  phi = phi + detGeo.getRingDPhi(iRing) / 2.;
1546  G4ThreeVector transBack(r * cos(phi), r * sin(phi), detGeo.getSupportThickness() / 2.);
1547  G4RotationMatrix RaBack;
1548  RaBack.rotateZ(phi);
1549  new G4PVPlacement(G4Transform3D(RaBack, transBack), hapdBackRadialWallLV[iRing - 1], "hapdBack", detSupportLV, false, iSlot);
1550 
1551  // add FEB cooling bodies
1552  G4ThreeVector febCoolingTa;
1553  G4Transform3D febCoolingTr;
1554  febCoolingTa.setX(r * cos(detGeo.getSlotPhi(iSlot)));
1555  febCoolingTa.setY(r * sin(detGeo.getSlotPhi(iSlot)));
1556 
1557  double supportTube_envelope_dZ = (detGeo.getSupportThickness() + detGeo.getSupportBackWallHeight());
1558  double febCooling_envelope_dZ = (detectorGeo.getFEBCoolingGeometry().getBigSquareThickness() +
1560  double febCooling_envelope_Z0 = -supportTube_envelope_dZ / 2.0 + febCooling_envelope_dZ / 2.0 + detGeo.getSupportThickness();
1561  febCoolingTa.setZ(febCooling_envelope_Z0);
1562 
1563  int febcoolingv2GeometryID = detectorGeo.getFEBCoolingGeometry().getFebcoolingv2GeometryID().at(iSlot - 1);
1564 
1565  if (febcoolingv2GeometryID == 2) Ra.rotateZ(90.0 * deg);
1566 
1567  febCoolingTr = G4Transform3D(Ra, febCoolingTa);
1568 
1569  if (febcoolingv2GeometryID != 0) {
1570  G4LogicalVolume* febCoolingLV = buildFEBCoolingBody(detectorGeo.getFEBCoolingGeometry());
1571 
1572  new G4PVPlacement(febCoolingTr, //Transformation
1573  febCoolingLV, //its logical volume
1574  "febCoolingLV", //its name
1575  detSupportLV, //its mother volume
1576  false, //no boolean operation
1577  iSlot); //copy number
1578  }
1579  }
1580 
1581  // G4LogicalVolume* hapdSupportPlateLV = new G4LogicalVolume(substraction, supportMaterial, "hapdSupport");
1582 
1583  G4Transform3D transform3 = G4Translate3D(0., 0., - backWallHeight / 2.);
1584  new G4PVPlacement(transform3, hapdSupportPlateLV, "supportPlate", detSupportLV, false, 1);
1585 
1586  // place electronics side neutron shield - for now hardcoded here, pending for more proper implementation!
1587  G4Box* shieldBox1 = new G4Box("shieldBox1", 20. / 2., 75. / 2., backWallHeight / 2.);
1588  G4Box* shieldBox2 = new G4Box("shieldBox2", 55. / 2., 40. / 2., backWallHeight / 2.);
1589  G4LogicalVolume* shield1 = new G4LogicalVolume(shieldBox1, Materials::get("BoratedPoly"), "ARICH.FWDShield1");
1590  G4LogicalVolume* shield2 = new G4LogicalVolume(shieldBox2, Materials::get("BoratedPoly"), "ARICH.FWDShield2");
1591  double dphi = 2 * M_PI / 36.;
1592  double r1 = wallR[0] - 15.;
1593  double r2 = wallR[0] - 15. - 20. / 2. - 55. / 2.;
1594  for (int i = 0; i < 36; i++) {
1595  double phi = (i + 0.5) * dphi;
1596  G4RotationMatrix rot;
1597  rot.rotateZ(phi);
1598  G4ThreeVector trans(r1 * cos(phi), r1 * sin(phi), detGeo.getSupportThickness() / 2.);
1599  G4ThreeVector trans1(r2 * cos(phi), r2 * sin(phi), detGeo.getSupportThickness() / 2.);
1600  new G4PVPlacement(G4Transform3D(rot, trans), shield1, "ARICH.FWDShield1", detSupportLV, false, i);
1601  new G4PVPlacement(G4Transform3D(rot, trans1), shield2, "ARICH.FWDShield2", detSupportLV, false, i);
1602  }
1603 
1604  return detSupportLV;
1605  }
1606 
1607 
1608  G4LogicalVolume* GeoARICHCreator::buildMirror(const ARICHGeometryConfig& detectorGeo)
1609  {
1610 
1611  const ARICHGeoMirrors& mirrGeo = detectorGeo.getMirrors();
1612 
1613  // read m_config
1614  string mirrMat = mirrGeo.getMaterial();
1615  G4Material* mirrorMaterial = Materials::get(mirrMat);
1616 
1617  double mThick = mirrGeo.getPlateThickness();
1618  double mLength = mirrGeo.getPlateLength();
1619  double mWidth = mirrGeo.getPlateWidth();
1620 
1621  G4Box* mirrPlate = new G4Box("mirrPlate", mThick / 2., mLength / 2., mWidth / 2.);
1622 
1623  G4LogicalVolume* lmirror = new G4LogicalVolume(mirrPlate, mirrorMaterial, "ARICH.mirrorPlate");
1624 
1625  Materials& materials = Materials::getInstance();
1626 
1627  G4OpticalSurface* optSurf = materials.createOpticalSurface(mirrGeo.getMirrorSurface());
1628  new G4LogicalSkinSurface("mirrorSurface", lmirror, optSurf);
1629 
1630  return lmirror;
1631 
1632  }
1633 
1634  double GeoARICHCreator::getAvgRINDEX(G4Material* material)
1635  {
1636  G4MaterialPropertiesTable* mTable = material->GetMaterialPropertiesTable();
1637  if (!mTable) return 0;
1638  G4MaterialPropertyVector* mVector = mTable->GetProperty("RINDEX");
1639  if (!mVector) return 0;
1640  G4bool b;
1641  return mVector->GetValue(2 * Unit::eV / Unit::MeV, b);
1642  }
1643 
1644  G4AssemblyVolume* GeoARICHCreator::makeJoint(G4Material* supportMaterial, const std::vector<double>& par)
1645  {
1646 
1647  int size = par.size();
1648  if (size < 4 || size > 8) B2ERROR("GeoARICHCreator::makeJoint: invalid number of joint wedge parameters");
1649  double lenx = par.at(0);
1650  double leny = par.at(1);
1651  double lenz = par.at(2);
1652  double thick = par.at(3);
1653 
1654  G4Box* wedgeBox1 = new G4Box("wedgeBox1", thick / 2., lenx / 2., leny / 2.);
1655  G4Box* wedgeBox2 = new G4Box("wedgeBox2", lenz / 2., lenx / 2., thick / 2.);
1656 
1657  G4LogicalVolume* wedgeBox1LV = new G4LogicalVolume(wedgeBox1, supportMaterial, "ARICH.supportWedge");
1658  G4LogicalVolume* wedgeBox2LV = new G4LogicalVolume(wedgeBox2, supportMaterial, "ARICH.supportWedge");
1659 
1660  G4AssemblyVolume* assemblyWedge = new G4AssemblyVolume();
1661 
1662  G4RotationMatrix Rm;
1663  G4ThreeVector Ta(0, 0, 0);
1664  G4Transform3D Tr;
1665  Tr = G4Transform3D(Rm, Ta);
1666 
1667  assemblyWedge->AddPlacedVolume(wedgeBox1LV, Tr);
1668 
1669  Ta.setX(lenz / 2. + thick / 2.);
1670  Ta.setZ(leny / 2. - thick / 2.);
1671  Tr = G4Transform3D(Rm, Ta);
1672  assemblyWedge->AddPlacedVolume(wedgeBox2LV, Tr);
1673 
1674  if (size == 4) return assemblyWedge;
1675  double edge = par.at(4);
1676 
1677  G4Box* wedgeBox3 = new G4Box("wedgeBox3", lenz / 2., edge / 2., thick / 2.);
1678  G4LogicalVolume* wedgeBox3LV = new G4LogicalVolume(wedgeBox3, supportMaterial, "ARICH.supportWedge");
1679 
1680  Ta.setZ(leny / 2. - thick - thick / 2.);
1681  Tr = G4Transform3D(Rm, Ta);
1682  assemblyWedge->AddPlacedVolume(wedgeBox3LV, Tr);
1683 
1684  G4Trap* wedgeBoxTmp = new G4Trap("wedgeBoxTmp", thick, leny - 2 * thick, lenz, edge);
1685  G4LogicalVolume* wedgeBox4LV;
1686  if (size == 8) {
1687  G4Tubs* wedgeBoxTmp1 = new G4Tubs("wedgeBoxTmp1", 0.0, par.at(5), thick, 0, 2.*M_PI);
1688  G4RotationMatrix rotHole;
1689  G4ThreeVector transHole(par.at(6), par.at(7), 0);
1690  G4SubtractionSolid* wedgeBox4 = new G4SubtractionSolid("wedgeBox4", wedgeBoxTmp, wedgeBoxTmp1, G4Transform3D(rotHole, transHole));
1691  wedgeBox4LV = new G4LogicalVolume(wedgeBox4, supportMaterial, "ARICH.supportWedge");
1692  } else wedgeBox4LV = new G4LogicalVolume(wedgeBoxTmp, supportMaterial, "ARICH.supportWedge");
1693 
1694  Rm.rotateX(-M_PI / 2.);
1695  Ta.setX(thick / 2. + edge / 4. + lenz / 4.);
1696  Ta.setZ(-thick / 2. - edge / 2.);
1697  Tr = G4Transform3D(Rm, Ta);
1698  assemblyWedge->AddPlacedVolume(wedgeBox4LV, Tr);
1699 
1700  return assemblyWedge;
1701 
1702  }
1703 
1704  // simple geometry for beamtest setup
1705  /* void GeoARICHCreator::createSimple(const GearDir& content, G4LogicalVolume& topVolume)
1706  {
1707 
1708  B2INFO("ARICH simple (beamtest) geometry will be built.")
1709  GearDir envelopeParams(content, "Envelope");
1710  double xSize = envelopeParams.getLength("xSize") / Unit::mm;
1711  double ySize = envelopeParams.getLength("ySize") / Unit::mm;
1712  double zSize = envelopeParams.getLength("zSize") / Unit::mm;
1713  string envMat = envelopeParams.getString("material");
1714  G4Material* envMaterial = Materials::get(envMat);
1715 
1716  G4Box* envBox = new G4Box("envBox", xSize / 2., ySize / 2., zSize / 2.);
1717  G4LogicalVolume* lenvBox = new G4LogicalVolume(envBox, envMaterial, "ARICH.envelope");
1718  new G4PVPlacement(G4Transform3D(), lenvBox, "ARICH.envelope", &topVolume, false, 1);
1719  setVisibility(*lenvBox, false);
1720  ARICHGeometryPar* m_arichgp = ARICHGeometryPar::Instance();
1721  m_arichgp->Initialize(content);
1722  GearDir moduleParam(content, "Detector/Module");
1723  G4LogicalVolume* detModule = buildModule(moduleParam);
1724 
1725  double detZpos = content.getLength("Detector/Plane/zPosition") / Unit::mm;
1726  double detThick = content.getLength("Detector/Module/moduleZSize") / Unit::mm;
1727  int nModules = m_arichgp->getNMCopies();
1728  for (int i = 1; i <= nModules; i++) {
1729  G4ThreeVector origin = m_arichgp->getOriginG4(i); origin.setZ(detZpos + detThick / 2.);
1730  double angle = m_arichgp->getModAngle(i);
1731  G4RotationMatrix Ra;
1732  Ra.rotateZ(angle);
1733  G4Transform3D trans = G4Transform3D(Ra, origin);
1734  new G4PVPlacement(G4Transform3D(Ra, origin), detModule, "detModule", lenvBox, false, i);
1735  }
1736 
1737  // place aerogel tiles
1738  GearDir aerogelParam(content, "Aerogel");
1739  double sizeX = aerogelParam.getLength("tileXSize") / Unit::mm;
1740  double sizeY = aerogelParam.getLength("tileYSize") / Unit::mm;
1741  double posX = aerogelParam.getLength("tileXPos") / Unit::mm;
1742  double posY = aerogelParam.getLength("tileYPos") / Unit::mm;
1743  int ilayer = 0;
1744  BOOST_FOREACH(const GearDir & layer, aerogelParam.getNodes("Layers/Layer")) {
1745  double posZ = layer.getLength("zPosition");
1746  double sizeZ = layer.getLength("thickness");
1747  string tileMat = layer.getString("material");
1748  G4Material* tileMaterial = Materials::get(tileMat);
1749  double rInd = getAvgRINDEX(tileMaterial);
1750  m_arichgp->setAeroRefIndex(ilayer, rInd);
1751  m_arichgp->setAerogelZPosition(ilayer, posZ);
1752  m_arichgp->setAerogelThickness(ilayer, sizeZ);
1753 
1754  G4MaterialPropertiesTable* mTable = tileMaterial->GetMaterialPropertiesTable();
1755  G4MaterialPropertyVector* mVector = mTable->GetProperty("RAYLEIGH");
1756  double lambda0 = 400;
1757  double e0 = 1240. / lambda0 * Unit::eV / Unit::MeV;
1758  G4bool b;
1759  m_arichgp->setAeroTransLength(ilayer, mVector->GetValue(e0, b)*Unit::mm);
1760  G4Box* tileBox = new G4Box("tileBox", sizeX / 2., sizeY / 2., sizeZ / 2. / Unit::mm);
1761  G4LogicalVolume* lTile = new G4LogicalVolume(tileBox, tileMaterial, "Tile", 0, m_sensitiveAero);
1762  setColor(*lTile, "rgb(0.0, 1.0, 1.0,1.0)");
1763  G4Transform3D trans = G4Translate3D(posX, posY, posZ / Unit::mm + sizeZ / 2. / Unit::mm);
1764  ilayer++;
1765  new G4PVPlacement(trans, lTile, "ARICH.tile", lenvBox, false, ilayer);
1766  }
1767 
1768  // place mirrors
1769  GearDir mirrorsParam(content, "Mirrors");
1770  double height = mirrorsParam.getLength("height") / Unit::mm;
1771  double width = mirrorsParam.getLength("width") / Unit::mm;
1772  double thickness = mirrorsParam.getLength("thickness") / Unit::mm;
1773  string mirrMat = mirrorsParam.getString("material");
1774  G4Material* mirrMaterial = Materials::get(mirrMat);
1775  G4Box* mirrBox = new G4Box("mirrBox", thickness / 2., height / 2., width / 2.);
1776  G4LogicalVolume* lmirror = new G4LogicalVolume(mirrBox, mirrMaterial, "mirror");
1777 
1778  Materials& materials = Materials::getInstance();
1779  GearDir surface(mirrorsParam, "Surface");
1780  G4OpticalSurface* optSurf = materials.createOpticalSurface(surface);
1781  new G4LogicalSkinSurface("mirrorsSurface", lmirror, optSurf);
1782  int iMirror = 0;
1783  BOOST_FOREACH(const GearDir & mirror, mirrorsParam.getNodes("Mirror")) {
1784  double xpos = mirror.getLength("xPos") / Unit::mm;
1785  double ypos = mirror.getLength("yPos") / Unit::mm;
1786  double zpos = mirror.getLength("zPos") / Unit::mm;
1787  double angle = mirror.getAngle("angle") / Unit::rad;
1788  G4ThreeVector origin(xpos, ypos, zpos + width / 2.);
1789  G4RotationMatrix Ra;
1790  Ra.rotateZ(angle);
1791  G4Transform3D trans = G4Transform3D(Ra, origin);
1792  new G4PVPlacement(G4Transform3D(Ra, origin), lmirror, "ARICH.mirror", lenvBox, false, iMirror);
1793  iMirror++;
1794  }
1795  }
1796  */
1797 
1798  } // namespace aricih
1800 } // 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
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
G4AssemblyVolume * makeJoint(G4Material *supportMaterial, const std::vector< double > &pars)
build joints of the ARICH support structure
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:286
int doBeamBackgroundStudy() const
returns 1 if beam background study (to add additional sensitive modules, detect neutrons,...
void setVisibility(G4LogicalVolume &volume, bool visible)
Helper function to quickly set the visibility of a given volume.
Definition: utilities.cc:105
void setColor(G4LogicalVolume &volume, const std::string &color)
Set the color of a logical volume.
Definition: utilities.cc:97
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.