Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 #include "HGCalTBMaterials.hh" 27 #include "G4Box.hh" 28 #include "G4Colour.hh" 29 #include "G4Element.hh" 30 #include "G4LogicalVolume.hh" 31 #include "G4Material.hh" 32 #include "G4NistManager.hh" 33 #include "G4PVPlacement.hh" 34 #include "G4RotationMatrix.hh" 35 #include "G4SubtractionSolid.hh" 36 #include "G4ThreeVector.hh" 37 #include "G4Tubs.hh" 38 #include "G4VisAttributes.hh" 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 G4SubtractionSolid *HexagonSolid(G4String aName, G4double aCellThickness, 43 G4double aCellSideLength) { 44 G4double fullCcellX = (2.) * aCellSideLength; 45 G4double fullCellY = std::sqrt(3.) * aCellSideLength; 46 G4Box *solidFullcell = new G4Box(aName, // its aName 47 0.5 * fullCcellX, 0.5 * fullCellY, 48 0.5 * aCellThickness); // its size 49 50 G4double deltaXDash = aCellSideLength; 51 G4double deltaYDash = std::sqrt(3) / 4 * aCellSideLength; 52 53 G4Box *solidCutcell = new G4Box(aName, // its aName 54 0.5 * deltaXDash, 0.5 * (deltaYDash), 55 1. * aCellThickness); // its size 56 57 G4double deltaTheta[4] = {30. * CLHEP::deg, 150. * CLHEP::deg, 210. * CLHEP::deg, 330. * CLHEP::deg}; 58 G4double deltaThetaRot[4] = {60. * CLHEP::deg, 120. * CLHEP::deg, 240 * CLHEP::deg, 300 * CLHEP::deg}; 59 G4double delta = std::sqrt(3) / 2 * aCellSideLength + deltaYDash / 2; 60 61 G4RotationMatrix *rot = new G4RotationMatrix; 62 rot->rotateZ(deltaThetaRot[0]); 63 std::vector<G4SubtractionSolid *> subtracted; 64 subtracted.push_back( 65 new G4SubtractionSolid("cellSubtracted", solidFullcell, solidCutcell, rot, 66 G4ThreeVector(std::cos(deltaTheta[0]) * delta, 67 std::sin(deltaTheta[0]) * delta, 0.))); 68 69 for (int i = 1; i < 4; i++) { 70 rot->rotateZ(-deltaThetaRot[i - 1]); 71 rot->rotateZ(deltaThetaRot[i]); 72 subtracted.push_back(new G4SubtractionSolid( 73 "cellSubtracted", subtracted[i - 1], solidCutcell, rot, 74 G4ThreeVector(std::cos(deltaTheta[i]) * delta, std::sin(deltaTheta[i]) * delta, 75 0.))); 76 } 77 delete rot; 78 79 return subtracted[3]; 80 } 81 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 83 84 G4LogicalVolume *HexagonLogical(G4String aName, G4double aCellThickness, 85 G4double aCellSideLength, 86 G4Material *aMaterial) { 87 return new G4LogicalVolume( 88 HexagonSolid(aName, aCellThickness, aCellSideLength), // its solid 89 aMaterial, // its aMaterial 90 aName); // its aName 91 } 92 93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 94 95 void HGCalTBMaterials::DefineMaterials() { 96 /***** Definition of all available materials *****/ 97 // Get nist aMaterial manager 98 G4NistManager *nist = G4NistManager::Instance(); 99 fMatVacuum = nist->FindOrBuildMaterial("G4_Galactic"); 100 fMatAIR = nist->FindOrBuildMaterial("G4_AIR"); 101 fMatAr = nist->FindOrBuildMaterial("G4_Ar"); 102 fMatAl = nist->FindOrBuildMaterial("G4_Al"); 103 fMatFe = nist->FindOrBuildMaterial("G4_Fe"); 104 fMatGlass = nist->FindOrBuildMaterial("G4_GLASS_PLATE"); 105 fMatPb = nist->FindOrBuildMaterial("G4_Pb"); 106 fMatCu = nist->FindOrBuildMaterial("G4_Cu"); 107 fMatW = nist->FindOrBuildMaterial("G4_W"); 108 fMatSi = nist->FindOrBuildMaterial("G4_Si"); 109 fMatAu = nist->FindOrBuildMaterial("G4_Au"); 110 fMatQuartz = nist->FindOrBuildMaterial("G4_SILICON_DIOXIDE"); 111 fMatC = nist->FindOrBuildMaterial("G4_C"); 112 fMatH = nist->FindOrBuildMaterial("G4_H"); 113 fMatO = nist->FindOrBuildMaterial("G4_O"); 114 fMatMn = nist->FindOrBuildMaterial("G4_Mn"); 115 fMatCr = nist->FindOrBuildMaterial("G4_Cr"); 116 fMatNi = nist->FindOrBuildMaterial("G4_Ni"); 117 fMatPolyethylene = nist->FindOrBuildMaterial("G4_POLYETHYLENE"); 118 fMatCl = nist->FindOrBuildMaterial("G4_Cl"); 119 fMatF = nist->FindOrBuildMaterial("G4_F"); 120 121 // AHCAL SiPMs 122 G4double a = 1.01 * CLHEP::g / CLHEP::mole; 123 G4Element *elH = new G4Element("Hydrogen", "H2", 1., a); 124 a = 12.01 * CLHEP::g / CLHEP::mole; 125 G4Element *elC = new G4Element("Carbon", "C", 6., a); 126 G4double density = 1.032 * CLHEP::g / CLHEP::cm3; 127 fMatPolystyrene = new G4Material("Polystyrene", density, 2); 128 fMatPolystyrene->AddElement(elC, 19); 129 fMatPolystyrene->AddElement(elH, 21); 130 131 // CuW alloy: 60% Cu, 40% W in mass 132 G4double CuFracInCuW = 0.75; 133 fMatCuW = new G4Material("CuW", 14.979 * CLHEP::g / CLHEP::cm3, 2); 134 fMatCuW->AddMaterial(fMatCu, CuFracInCuW); 135 fMatCuW->AddMaterial(fMatW, 1 - CuFracInCuW); 136 137 // PCB aMaterial 138 fMatPCB = new G4Material("PCB", 1.7 * CLHEP::g / CLHEP::cm3, 4); 139 fMatPCB->AddMaterial(fMatC, 0.13232243); 140 fMatPCB->AddMaterial(fMatH, 0.032572448); 141 fMatPCB->AddMaterial(fMatO, 0.48316123); 142 fMatPCB->AddMaterial(fMatSi, 0.35194389); 143 144 // Kapton aMaterial 145 fMatKAPTON = new G4Material("Kapton", 1.11 * CLHEP::g / CLHEP::cm3, 3); 146 fMatKAPTON->AddMaterial(fMatC, 0.59985105); 147 fMatKAPTON->AddMaterial(fMatH, 0.080541353); 148 fMatKAPTON->AddMaterial(fMatO, 0.31960759); 149 150 // steel 151 fMatSteel = new G4Material("StainlessSteel", 8.02 * CLHEP::g / CLHEP::cm3, 5); 152 fMatSteel->AddMaterial(fMatFe, 0.6996); 153 fMatSteel->AddMaterial(fMatC, 0.0004); 154 fMatSteel->AddMaterial(fMatMn, 0.01); 155 fMatSteel->AddMaterial(fMatCr, 0.19); 156 fMatSteel->AddMaterial(fMatNi, 0.1); 157 158 // Scintillator aMaterial 159 fMatScintillator = new G4Material("Scintillator", 1.032 * CLHEP::g / CLHEP::cm3, 2); 160 fMatScintillator->AddMaterial(fMatC, 0.91512109); 161 fMatScintillator->AddMaterial(fMatH, 0.084878906); 162 163 // DWC gas 164 fMatArCO2 = new G4Material("ArCO2", 1.729 * CLHEP::mg / CLHEP::cm3, 3); 165 fMatArCO2->AddMaterial(fMatAr, 0.475815); 166 fMatArCO2->AddMaterial(fMatC, 0.14306133); 167 fMatArCO2->AddMaterial(fMatO, 0.38112367); 168 169 fMatFreon = new G4Material("Freon-12", 4.93 * CLHEP::mg / CLHEP::cm3, 3); 170 fMatFreon->AddMaterial(fMatC, 0.099340816); 171 fMatFreon->AddMaterial(fMatCl, 0.58640112); 172 fMatFreon->AddMaterial(fMatF, 0.31425807); 173 } 174 175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 176 177 void HGCalTBMaterials::SetEventDisplayColorScheme() { 178 G4VisAttributes *visAttributes; 179 180 visAttributes = new G4VisAttributes(G4Colour(.0, 0.0, 0.0)); 181 visAttributes->SetVisibility(false); 182 fSiWaferLogical->SetVisAttributes(visAttributes); 183 184 visAttributes = new G4VisAttributes(G4Colour(.3, 0.3, 0.3, 0.2)); 185 visAttributes->SetVisibility(true); 186 fSiPixelLogical->SetVisAttributes(visAttributes); 187 188 visAttributes = new G4VisAttributes(G4Colour(0.4, 0.4, 0.4, 0.01)); 189 visAttributes->SetVisibility(false); 190 fCuWbaseplateLogical->SetVisAttributes(visAttributes); 191 fCuWbaseplate550umLogical->SetVisAttributes(visAttributes); 192 fCuWbaseplate610umLogical->SetVisAttributes(visAttributes); 193 fCuWbaseplate710umLogical->SetVisAttributes(visAttributes); 194 fCuBaseplateLogical->SetVisAttributes(visAttributes); 195 fCuBaseplate25umLogical->SetVisAttributes(visAttributes); 196 fCuBaseplate175umLogical->SetVisAttributes(visAttributes); 197 fPCBbaseplateLogical->SetVisAttributes(visAttributes); 198 fPCBbaseplateThinLogical->SetVisAttributes(visAttributes); 199 fKaptonLayerLogical->SetVisAttributes(visAttributes); 200 fAlCaseLogical->SetVisAttributes(visAttributes); 201 fAlCaseThickLogical->SetVisAttributes(visAttributes); 202 fSteelCaseLogical->SetVisAttributes(visAttributes); 203 fSteelCaseThickLogical->SetVisAttributes(visAttributes); 204 fPbAbsorberEElogical->SetVisAttributes(visAttributes); 205 fFeAbsorberEElogical->SetVisAttributes(visAttributes); 206 fCuAbsorberEElogical->SetVisAttributes(visAttributes); 207 fWabsorberEElogical->SetVisAttributes(visAttributes); 208 fW2mmAbsorberEEDESY2018Logical->SetVisAttributes(visAttributes); 209 fW4mmAbsorberEEDESY2018Logical->SetVisAttributes(visAttributes); 210 fCuAbsorberFHlogical->SetVisAttributes(visAttributes); 211 fFeAbsorberFHlogical->SetVisAttributes(visAttributes); 212 fAHCALSiPMlogical->SetVisAttributes(visAttributes); 213 fAHCALSiPM2x2HUBlogical->SetVisAttributes(visAttributes); 214 fAlAbsorberAHCALlogical->SetVisAttributes(visAttributes); 215 fPCBAHCALlogical->SetVisAttributes(visAttributes); 216 fFeAbsorberAHCALlogical->SetVisAttributes(visAttributes); 217 fScintillatorLogical->SetVisAttributes(visAttributes); 218 fScintillatorThinLogical->SetVisAttributes(visAttributes); 219 fMCPlogical->SetVisAttributes(visAttributes); 220 fCK3logical->SetVisAttributes(visAttributes); 221 fDWClogical->SetVisAttributes(visAttributes); 222 fDWCgasLogical->SetVisAttributes(visAttributes); 223 fAlChipLogical->SetVisAttributes(visAttributes); 224 } 225 226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 227 228 HGCalTBMaterials::HGCalTBMaterials() { 229 DefineMaterials(); 230 DefineSiWaferAndCells(); 231 DefineHGCalBaseplates(); 232 DefineHGCalCases(); 233 DefineHGCalEEAbsorbers(); 234 DefineHGCalFHAbsorbers(); 235 DefineAHCALSiPM(); 236 DefineAHCALAbsorbers(); 237 DefineBeamLineElements(); 238 } 239 240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 241 242 void HGCalTBMaterials::DefineSiWaferAndCells() { 243 /***** Definition of silicon (wafer) sensors *****/ 244 // 300 microns thickness only 245 fSiPixelSideLength = 0.6496345 * CLHEP::cm; 246 fSiWaferThickness = 0.3 * CLHEP::mm; 247 fAlpha = 60. / 180. * CLHEP::pi; 248 fSiWaferSideLength = 11 * fSiPixelSideLength; 249 250 fSiWaferLogical = HexagonLogical("Si_wafer", fSiWaferThickness, 251 fSiWaferSideLength, fMatAIR); 252 253 // Silicon pixel setups 254 double dx = 2 * std::sin(fAlpha) * fSiPixelSideLength; 255 double dy = fSiPixelSideLength * (2. + 2 * std::cos(fAlpha)); 256 fSiPixelLogical = 257 HexagonLogical("SiCell", fSiWaferThickness, fSiPixelSideLength, fMatSi); 258 259 int index = 1; 260 int nRows[11] = {7, 6, 7, 6, 5, 6, 5, 4, 5, 4, 3}; 261 for (int nC = 0; nC < 11; nC++) { 262 for (int middle_index = 0; middle_index < nRows[nC]; middle_index++) { 263 new G4PVPlacement( 264 0, 265 G4ThreeVector(dy * (middle_index - nRows[nC] / 2. + 0.5), nC * dx / 2, 266 0.), 267 fSiPixelLogical, "SiCell", fSiWaferLogical, false, index++, true); 268 if (nC <= 0) 269 continue; 270 new G4PVPlacement( 271 0, 272 G4ThreeVector(dy * (middle_index - nRows[nC] / 2. + 0.5), 273 -nC * dx / 2, 0.), 274 fSiPixelLogical, "SiCell", fSiWaferLogical, false, index++, true); 275 } 276 } 277 278 fThicknessMap["Si_wafer"] = fSiWaferThickness; 279 fLogicalVolumeMap["Si_wafer"] = fSiWaferLogical; 280 } 281 282 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 283 284 void HGCalTBMaterials::DefineHGCalBaseplates() { 285 /***** Definition of all baseplates *****/ 286 // CuW 287 G4double CuWbaseplateThickness = 1.2 * CLHEP::mm; 288 G4double CuWbaseplatesideLength = 11 * fSiPixelSideLength; 289 fCuWbaseplateLogical = HexagonLogical("CuW_baseplate", CuWbaseplateThickness, 290 CuWbaseplatesideLength, fMatCuW); 291 fThicknessMap["CuW_baseplate"] = CuWbaseplateThickness; 292 fLogicalVolumeMap["CuW_baseplate"] = fCuWbaseplateLogical; 293 G4double CuWbaseplate550umthickness = 0.55 * CLHEP::mm; 294 fCuWbaseplate550umLogical = 295 HexagonLogical("CuW_baseplate_550um", CuWbaseplate550umthickness, 296 CuWbaseplatesideLength, fMatCuW); 297 fThicknessMap["CuW_baseplate_550um"] = CuWbaseplate550umthickness; 298 fLogicalVolumeMap["CuW_baseplate_550um"] = fCuWbaseplate550umLogical; 299 G4double CuWbaseplate610umThickness = 0.61 * CLHEP::mm; 300 fCuWbaseplate610umLogical = 301 HexagonLogical("CuW_baseplate_610um", CuWbaseplate610umThickness, 302 CuWbaseplatesideLength, fMatCuW); 303 fThicknessMap["CuW_baseplate_610um"] = CuWbaseplate610umThickness; 304 fLogicalVolumeMap["CuW_baseplate_610um"] = fCuWbaseplate610umLogical; 305 G4double CuWbaseplate710umThickness = 0.71 * CLHEP::mm; 306 fCuWbaseplate710umLogical = 307 HexagonLogical("CuW_baseplate_710um", CuWbaseplate710umThickness, 308 CuWbaseplatesideLength, fMatCuW); 309 fThicknessMap["CuW_baseplate_710um"] = CuWbaseplate710umThickness; 310 fLogicalVolumeMap["CuW_baseplate_710um"] = fCuWbaseplate710umLogical; 311 312 // Cu 313 G4double CuBaseplateThickness = 1.2 * CLHEP::mm; 314 G4double CuBaseplatesideLength = 11 * fSiPixelSideLength; 315 fCuBaseplateLogical = HexagonLogical("Cu_baseplate", CuBaseplateThickness, 316 CuBaseplatesideLength, fMatCu); 317 fThicknessMap["Cu_baseplate"] = CuBaseplateThickness; 318 fLogicalVolumeMap["Cu_baseplate"] = fCuBaseplateLogical; 319 G4double CuBaseplate25umThickness = 0.025 * CLHEP::mm; 320 fCuBaseplate25umLogical = 321 HexagonLogical("Cu_baseplate_25um", CuBaseplate25umThickness, 322 CuBaseplatesideLength, fMatCu); 323 fThicknessMap["Cu_baseplate_25um"] = CuBaseplate25umThickness; 324 fLogicalVolumeMap["Cu_baseplate_25um"] = fCuBaseplate25umLogical; 325 G4double CuBaseplate175umThickness = 0.175 * CLHEP::mm; 326 fCuBaseplate175umLogical = 327 HexagonLogical("Cu_baseplate_175um", CuBaseplate175umThickness, 328 CuBaseplatesideLength, fMatCu); 329 fThicknessMap["Cu_baseplate_175um"] = CuBaseplate175umThickness; 330 fLogicalVolumeMap["Cu_baseplate_175um"] = fCuBaseplate175umLogical; 331 332 // PCB 333 G4double PCBbaseplateThickness = 1.3 * CLHEP::mm; 334 G4double PCBbaseplateSideLength = 11 * fSiPixelSideLength; 335 fPCBbaseplateLogical = HexagonLogical("PCB", PCBbaseplateThickness, 336 PCBbaseplateSideLength, fMatPCB); 337 fThicknessMap["PCB"] = PCBbaseplateThickness; 338 fLogicalVolumeMap["PCB"] = fPCBbaseplateLogical; 339 340 G4double PCBbaseplateThinThickness = 1.2 * CLHEP::mm; 341 fPCBbaseplateThinLogical = HexagonLogical( 342 "PCB_thin", PCBbaseplateThinThickness, PCBbaseplateSideLength, fMatPCB); 343 fThicknessMap["PCB_thin"] = PCBbaseplateThinThickness; 344 fLogicalVolumeMap["PCB_thin"] = fPCBbaseplateThinLogical; 345 346 // Kapton layer 347 G4double KaptonLayerThickness = 0.075 * CLHEP::mm; 348 G4double KaptonLayerSideLength = 11 * fSiPixelSideLength; 349 fKaptonLayerLogical = HexagonLogical("Kapton_layer", KaptonLayerThickness, 350 KaptonLayerSideLength, fMatKAPTON); 351 fThicknessMap["Kapton_layer"] = KaptonLayerThickness; 352 fLogicalVolumeMap["Kapton_layer"] = fKaptonLayerLogical; 353 } 354 355 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 356 357 void HGCalTBMaterials::DefineHGCalCases() { 358 G4double AlCaseThickness = 2.1 * CLHEP::mm; 359 G4double AlCaseXY = 40 * CLHEP::cm; 360 G4Box *AlCaseSolid = new G4Box("Al_case", 0.5 * AlCaseXY, 0.5 * AlCaseXY, 361 0.5 * AlCaseThickness); 362 fAlCaseLogical = new G4LogicalVolume(AlCaseSolid, fMatAl, "Al_case"); 363 fThicknessMap["Al_case"] = AlCaseThickness; 364 fLogicalVolumeMap["Al_case"] = fAlCaseLogical; 365 366 G4double AlCaseThickThickness = 5 * CLHEP::mm; 367 G4Box *AlCaseThickSolid = 368 new G4Box("Al_case_thick", 0.5 * AlCaseXY, 0.5 * AlCaseXY, 369 0.5 * AlCaseThickThickness); 370 fAlCaseThickLogical = 371 new G4LogicalVolume(AlCaseThickSolid, fMatAl, "Al_case_thick"); 372 fThicknessMap["Al_case_thick"] = AlCaseThickThickness; 373 fLogicalVolumeMap["Al_case_thick"] = fAlCaseThickLogical; 374 375 G4double SteelCaseThickness = 9 * CLHEP::mm; 376 G4double SteelCaseXY = 60 * CLHEP::cm; 377 G4Box *SteelCaseSolid = 378 new G4Box("Steel_case", 0.5 * SteelCaseXY, 0.5 * SteelCaseXY, 379 0.5 * SteelCaseThickness); 380 fSteelCaseLogical = new G4LogicalVolume(SteelCaseSolid, fMatFe, "Steel_case"); 381 fThicknessMap["Steel_case"] = SteelCaseThickness; 382 fLogicalVolumeMap["Steel_case"] = fSteelCaseLogical; 383 384 G4double SteelCaseThickThickness = 40 * CLHEP::mm; 385 G4Box *SteelCaseThickSolid = 386 new G4Box("Steel_case_thick", 0.5 * SteelCaseXY, 0.5 * SteelCaseXY, 387 0.5 * SteelCaseThickThickness); 388 fSteelCaseThickLogical = 389 new G4LogicalVolume(SteelCaseThickSolid, fMatFe, "Steel_case_thick"); 390 fThicknessMap["Steel_case_thick"] = SteelCaseThickThickness; 391 fLogicalVolumeMap["Steel_case_thick"] = fSteelCaseThickLogical; 392 } 393 394 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 395 396 void HGCalTBMaterials::DefineHGCalEEAbsorbers() { 397 // defintion of absorber plates in the EE part 398 G4double PbAbsorberEEthickness = 4.9 * CLHEP::mm; 399 G4double PbAbsorberEExy = 30 * CLHEP::cm; 400 G4Box *PbAbsorberEEsolid = 401 new G4Box("Pb_absorber_EE", 0.5 * PbAbsorberEExy, 0.5 * PbAbsorberEExy, 402 0.5 * PbAbsorberEEthickness); 403 fPbAbsorberEElogical = 404 new G4LogicalVolume(PbAbsorberEEsolid, fMatPb, "Pb_absorber_EE"); 405 fThicknessMap["Pb_absorber_EE"] = PbAbsorberEEthickness; 406 fLogicalVolumeMap["Pb_absorber_EE"] = fPbAbsorberEElogical; 407 408 G4double FeAbsorberEEthickness = 0.3 * CLHEP::mm; 409 G4double FeAbsorberEExy = 30 * CLHEP::cm; 410 G4Box *FeAbsorberEEsolid = 411 new G4Box("Fe_absorber_EE", 0.5 * FeAbsorberEExy, 0.5 * FeAbsorberEExy, 412 0.5 * FeAbsorberEEthickness); 413 fFeAbsorberEElogical = 414 new G4LogicalVolume(FeAbsorberEEsolid, fMatFe, "Fe_absorber_EE"); 415 fThicknessMap["Fe_absorber_EE"] = FeAbsorberEEthickness; 416 fLogicalVolumeMap["Fe_absorber_EE"] = fFeAbsorberEElogical; 417 418 G4double CuAbsorberEEthickness = 6 * CLHEP::mm; 419 G4double CuAbsorberEExy = 30 * CLHEP::cm; 420 G4Box *CuAbsorberEEsolid = 421 new G4Box("Cu_absorber_EE", 0.5 * CuAbsorberEExy, 0.5 * CuAbsorberEExy, 422 0.5 * CuAbsorberEEthickness); 423 fCuAbsorberEElogical = 424 new G4LogicalVolume(CuAbsorberEEsolid, fMatCu, "Cu_absorber_EE"); 425 fThicknessMap["Cu_absorber_EE"] = CuAbsorberEEthickness; 426 fLogicalVolumeMap["Cu_absorber_EE"] = fCuAbsorberEElogical; 427 428 G4double WabsorberEEthickness = 2.8 * CLHEP::mm; 429 G4double WabsorberEExy = 30 * CLHEP::cm; 430 G4Box *WabsorberEEsolid = 431 new G4Box("W_absorber_EE", 0.5 * WabsorberEExy, 0.5 * WabsorberEExy, 432 0.5 * WabsorberEEthickness); 433 fWabsorberEElogical = 434 new G4LogicalVolume(WabsorberEEsolid, fMatW, "W_absorber_EE"); 435 fThicknessMap["W_absorber_EE"] = WabsorberEEthickness; 436 fLogicalVolumeMap["W_absorber_EE"] = fWabsorberEElogical; 437 438 G4double W2mmAbsorberEEDESY2018thickness = 2. * CLHEP::mm; 439 G4double W2mmAbsorberEEDESY2018xy = 15 * CLHEP::cm; 440 G4Box *W2mmAbsorberEEDESY2018solid = new G4Box( 441 "W_2mm_absorber_EE_DESY2018", 0.5 * W2mmAbsorberEEDESY2018xy, 442 0.5 * W2mmAbsorberEEDESY2018xy, 0.5 * W2mmAbsorberEEDESY2018thickness); 443 fW2mmAbsorberEEDESY2018Logical = new G4LogicalVolume( 444 W2mmAbsorberEEDESY2018solid, fMatW, "W_2mm_absorber_EE_DESY2018"); 445 fThicknessMap["W_2mm_absorber_EE_DESY2018"] = W2mmAbsorberEEDESY2018thickness; 446 fLogicalVolumeMap["W_2mm_absorber_EE_DESY2018"] = 447 fW2mmAbsorberEEDESY2018Logical; 448 449 G4double W4mmAbsorberEEDESY2018thickness = 4. * CLHEP::mm; 450 G4double W4mmAbsorberEEDESY2018xy = 15 * CLHEP::cm; 451 G4Box *W4mmAbsorberEEDESY2018solid = new G4Box( 452 "W_4mm_absorber_EE_DESY2018", 0.5 * W4mmAbsorberEEDESY2018xy, 453 0.5 * W4mmAbsorberEEDESY2018xy, 0.5 * W4mmAbsorberEEDESY2018thickness); 454 fW4mmAbsorberEEDESY2018Logical = new G4LogicalVolume( 455 W4mmAbsorberEEDESY2018solid, fMatW, "W_4mm_absorber_EE_DESY2018"); 456 fThicknessMap["W_4mm_absorber_EE_DESY2018"] = W4mmAbsorberEEDESY2018thickness; 457 fLogicalVolumeMap["W_4mm_absorber_EE_DESY2018"] = 458 fW4mmAbsorberEEDESY2018Logical; 459 } 460 461 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 462 463 void HGCalTBMaterials::DefineHGCalFHAbsorbers() { 464 // defintion of absorber plates in the FH part 465 G4double CuAbsorberFHthickness = 6 * CLHEP::mm; 466 G4double CuAbsorberFHxy = 50 * CLHEP::cm; 467 G4Box *CuAbsorberFHsolid = 468 new G4Box("Cu_absorber_FH", 0.5 * CuAbsorberFHxy, 0.5 * CuAbsorberFHxy, 469 0.5 * CuAbsorberFHthickness); 470 fCuAbsorberFHlogical = 471 new G4LogicalVolume(CuAbsorberFHsolid, fMatCu, "Cu_absorber_FH"); 472 fThicknessMap["Cu_absorber_FH"] = CuAbsorberFHthickness; 473 fLogicalVolumeMap["Cu_absorber_FH"] = fCuAbsorberFHlogical; 474 475 G4double FeAbsorberFHthickness = 40 * CLHEP::mm; 476 G4double FeAbsorberFHxy = 50 * CLHEP::cm; 477 G4Box *FeAbsorberFHsolid = 478 new G4Box("Fe_absorber_FH", 0.5 * FeAbsorberFHxy, 0.5 * FeAbsorberFHxy, 479 0.5 * FeAbsorberFHthickness); 480 fFeAbsorberFHlogical = 481 new G4LogicalVolume(FeAbsorberFHsolid, fMatFe, "Fe_absorber_FH"); 482 fThicknessMap["Fe_absorber_FH"] = FeAbsorberFHthickness; 483 fLogicalVolumeMap["Fe_absorber_FH"] = fFeAbsorberFHlogical; 484 } 485 486 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 487 488 void HGCalTBMaterials::DefineAHCALSiPM() { 489 G4double AHCALSiPMthickness = 5.4 * CLHEP::mm; 490 fAHCALSiPMxy = 3 * CLHEP::cm; 491 fAHCALSiPMsolid = new G4Box("AHCAL_SiPM", 0.5 * fAHCALSiPMxy, 492 0.5 * fAHCALSiPMxy, 0.5 * AHCALSiPMthickness); 493 fAHCALSiPMlogical = 494 new G4LogicalVolume(fAHCALSiPMsolid, fMatPolystyrene, "AHCAL_SiPM"); 495 fThicknessMap["AHCAL_SiPM"] = AHCALSiPMthickness; 496 fLogicalVolumeMap["AHCAL_SiPM"] = fAHCALSiPMlogical; 497 498 G4double AHCALSiPM2x2HUBxy = 2 * 12 * fAHCALSiPMxy + 0.01 * CLHEP::mm; 499 G4double AHCALSiPM2x2HUBthickness = AHCALSiPMthickness + 0.01 * CLHEP::mm; 500 G4Box *AHCALSiPM2x2HUBsolid = 501 new G4Box("AHCAL_SiPM_2x2HUB", 0.5 * AHCALSiPM2x2HUBxy, 502 0.5 * AHCALSiPM2x2HUBxy, 0.5 * AHCALSiPM2x2HUBthickness); 503 fAHCALSiPM2x2HUBlogical = 504 new G4LogicalVolume(AHCALSiPM2x2HUBsolid, fMatAIR, "AHCAL_SiPM_2x2HUB"); 505 fThicknessMap["AHCAL_SiPM_2x2HUB"] = AHCALSiPM2x2HUBthickness; 506 fLogicalVolumeMap["AHCAL_SiPM_2x2HUB"] = fAHCALSiPM2x2HUBlogical; 507 int copy_counter = 0; 508 for (float _dx = -11.5; _dx <= 11.5; _dx = _dx + 1.) 509 for (float _dy = -11.5; _dy <= 11.5; _dy = _dy + 1.) 510 new G4PVPlacement( 511 0, G4ThreeVector(_dx * fAHCALSiPMxy, _dy * fAHCALSiPMxy, 0), 512 fAHCALSiPMlogical, "AHCAL_SiPM", fAHCALSiPM2x2HUBlogical, false, 513 copy_counter++, true); 514 } 515 516 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 517 518 void HGCalTBMaterials::DefineAHCALAbsorbers() { 519 G4double AlAbsorberAHCALthickness = 1 * CLHEP::mm; 520 G4double AlAbsorberAHCALxy = 2 * 12 * fAHCALSiPMxy; 521 G4Box *AlAbsorberAHCALsolid = 522 new G4Box("Al_absorber_AHCAL", 0.5 * AlAbsorberAHCALxy, 523 0.5 * AlAbsorberAHCALxy, 0.5 * AlAbsorberAHCALthickness); 524 fAlAbsorberAHCALlogical = 525 new G4LogicalVolume(AlAbsorberAHCALsolid, fMatAl, "Al_absorber_AHCAL"); 526 fThicknessMap["Al_absorber_AHCAL"] = AlAbsorberAHCALthickness; 527 fLogicalVolumeMap["Al_absorber_AHCAL"] = fAlAbsorberAHCALlogical; 528 529 G4double PCBAHCALthickness = 1.2 * CLHEP::mm; 530 G4double PCBAHCALxy = 2 * 12 * fAHCALSiPMxy; 531 G4Box *PCBAHCALsolid = new G4Box("PCB_AHCAL", 0.5 * PCBAHCALxy, 532 0.5 * PCBAHCALxy, 0.5 * PCBAHCALthickness); 533 fPCBAHCALlogical = new G4LogicalVolume(PCBAHCALsolid, fMatPCB, "PCB_AHCAL"); 534 fThicknessMap["PCB_AHCAL"] = PCBAHCALthickness; 535 fLogicalVolumeMap["PCB_AHCAL"] = fPCBAHCALlogical; 536 537 G4double FeAbsorberAHCALthickness = 17 * CLHEP::mm; 538 G4double FeAbsorberAHCALx = 80.8 * CLHEP::cm; 539 G4double FeAbsorberAHCALy = 65.7 * CLHEP::cm; 540 G4Box *FeAbsorberAHCALsolid = 541 new G4Box("Fe_absorber_AHCAL", 0.5 * FeAbsorberAHCALx, 542 0.5 * FeAbsorberAHCALy, 0.5 * FeAbsorberAHCALthickness); 543 fFeAbsorberAHCALlogical = 544 new G4LogicalVolume(FeAbsorberAHCALsolid, fMatFe, "Fe_absorber_AHCAL"); 545 fThicknessMap["Fe_absorber_AHCAL"] = FeAbsorberAHCALthickness; 546 fLogicalVolumeMap["Fe_absorber_AHCAL"] = fFeAbsorberAHCALlogical; 547 } 548 549 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 550 551 void HGCalTBMaterials::DefineBeamLineElements() { 552 /***** Definition of beam line elements *****/ 553 // scintillators 554 G4double scintillatorThickness = 2 * 2 * CLHEP::cm; 555 G4double scintillatorXY = 2 * 9.5 * CLHEP::cm; 556 G4Box *scintillatorSolid = 557 new G4Box("Scintillator", 0.5 * scintillatorXY, 0.5 * scintillatorXY, 558 0.5 * scintillatorThickness); 559 fScintillatorLogical = 560 new G4LogicalVolume(scintillatorSolid, fMatScintillator, "Scintillator"); 561 fThicknessMap["Scintillator"] = scintillatorThickness; 562 fLogicalVolumeMap["Scintillator"] = fScintillatorLogical; 563 564 G4double fScintillatorThinThickness = 1 * CLHEP::cm; 565 G4Box *fScintillatorThinSolid = 566 new G4Box("Scintillator_thin", 0.5 * scintillatorXY, 0.5 * scintillatorXY, 567 0.5 * fScintillatorThinThickness); 568 fScintillatorThinLogical = new G4LogicalVolume( 569 fScintillatorThinSolid, fMatScintillator, "Scintillator_thin"); 570 fThicknessMap["Scintillator_thin"] = fScintillatorThinThickness; 571 fLogicalVolumeMap["Scintillator_thin"] = fScintillatorThinLogical; 572 573 // MCPs = quartz disks 574 G4double MCPthickness = 10 * CLHEP::mm; 575 G4double MCPradius = 2 * CLHEP::cm; 576 G4Tubs *MCPsolid = 577 new G4Tubs("MCP", 0., MCPradius, MCPthickness, 0, 360 * CLHEP::degree); 578 fMCPlogical = new G4LogicalVolume(MCPsolid, fMatQuartz, "MCP"); 579 fThicknessMap["MCP"] = MCPthickness; 580 fLogicalVolumeMap["MCP"] = fMCPlogical; 581 582 // CK3 583 G4double CK3thickness = 2 * CLHEP::m; 584 G4double CK3radius = 8.35 * CLHEP::cm; 585 G4Tubs *CK3solid = 586 new G4Tubs("CK3", 0., CK3radius, 0.5 * CK3thickness, 0, 360 * CLHEP::degree); 587 fCK3logical = new G4LogicalVolume(CK3solid, fMatFreon, "CK3"); 588 fThicknessMap["CK3"] = CK3thickness; 589 fLogicalVolumeMap["CK3"] = fCK3logical; 590 591 // Aluminium circle for testing of chip impact 592 G4double AlChipXY = 1 * CLHEP::cm; 593 G4double AlChipThickness = 0.89 * CLHEP::mm; // corresponds to 1% X0 594 G4Box *AlChipSolid = new G4Box("Al_chip", 0.5 * AlChipXY, 0.5 * AlChipXY, 595 0.5 * AlChipThickness); 596 fAlChipLogical = new G4LogicalVolume(AlChipSolid, fMatAl, "Al_chip"); 597 fThicknessMap["Al_chip"] = AlChipThickness; 598 fLogicalVolumeMap["Al_chip"] = fAlChipLogical; 599 600 // DWC related aMaterial 601 G4double DWCthickness = 2 * 27.5 * CLHEP::mm; 602 G4double DWCxy = 2 * 11 * CLHEP::cm; 603 G4Box *DWCsolid = 604 new G4Box("DWC", 0.5 * DWCxy, 0.5 * DWCxy, 0.5 * DWCthickness); 605 fDWClogical = new G4LogicalVolume(DWCsolid, fMatAIR, "DWC"); 606 fThicknessMap["DWC"] = DWCthickness; 607 fLogicalVolumeMap["DWC"] = fDWClogical; 608 609 // WChambGas 610 G4double DWCgasThickness = 2 * 22.5 * CLHEP::mm; 611 G4double DWCgasXY = 2 * 8.5 * CLHEP::cm; 612 G4Box *DWCgasSolid = new G4Box("DWC_gas", 0.5 * DWCgasXY, 0.5 * DWCgasXY, 613 0.5 * DWCgasThickness); 614 fDWCgasLogical = new G4LogicalVolume(DWCgasSolid, fMatArCO2, "DWC_gas"); 615 new G4PVPlacement(0, G4ThreeVector(0, 0., 0.), fDWCgasLogical, "DWC_gas", 616 fDWClogical, false, 0, true); 617 618 G4VisAttributes *visAttributes = new G4VisAttributes(G4Colour(.0, 0.0, 0.0)); 619 visAttributes->SetVisibility(false); 620 // WChambWindow 621 G4double DWCwindowThickness = 0.025 * CLHEP::mm; 622 G4double DWCwindowXY = 2 * 5.5 * CLHEP::cm; 623 G4Box *DWCwindowSolid = 624 new G4Box("DWC_window", 0.5 * DWCwindowXY, 0.5 * DWCwindowXY, 625 0.5 * DWCwindowThickness); 626 auto DWCwindowLogical = 627 new G4LogicalVolume(DWCwindowSolid, fMatKAPTON, "DWC_window"); 628 DWCwindowLogical->SetVisAttributes(visAttributes); 629 new G4PVPlacement( 630 0, G4ThreeVector(0, 0., 27.5 * CLHEP::mm - DWCwindowThickness / 2.), 631 DWCwindowLogical, "DWC_window_0", fDWClogical, true, 0, true); 632 new G4PVPlacement( 633 0, G4ThreeVector(0, 0., -(27.5 * CLHEP::mm - DWCwindowThickness / 2.)), 634 DWCwindowLogical, "DWC_window_1", fDWClogical, true, 1, true); 635 636 G4RotationMatrix *rotation = new G4RotationMatrix(); 637 rotation->rotateZ(90 * CLHEP::deg); 638 639 // WChambAl1 640 G4double DWCal1thickness = 2 * 2.5 * CLHEP::mm; 641 G4double DWCal1x = 2 * 8.25 * CLHEP::cm; 642 G4double DWCal1y = 2 * 2.75 * CLHEP::cm; 643 G4Box *DWCal1solid = 644 new G4Box("DWC_al1", 0.5 * DWCal1x, 0.5 * DWCal1y, 0.5 * DWCal1thickness); 645 auto DWCal1logical = new G4LogicalVolume(DWCal1solid, fMatAl, "DWC_al1"); 646 DWCal1logical->SetVisAttributes(visAttributes); 647 new G4PVPlacement(0, 648 G4ThreeVector(0.5 * DWCal1y, 0.5 * DWCal1x, 649 (0.5 * DWCthickness - 0.5 * DWCal1thickness)), 650 DWCal1logical, "DWC_al1_0", fDWClogical, true, 0, true); 651 new G4PVPlacement(rotation, 652 G4ThreeVector(-0.5 * DWCal1x, 0.5 * DWCal1y, 653 (0.5 * DWCthickness - 0.5 * DWCal1thickness)), 654 DWCal1logical, "DWC_al1_1", fDWClogical, true, 1, true); 655 new G4PVPlacement(0, 656 G4ThreeVector(-0.5 * DWCal1y, -0.5 * DWCal1x, 657 (0.5 * DWCthickness - 0.5 * DWCal1thickness)), 658 DWCal1logical, "DWC_al1_2", fDWClogical, true, 2, true); 659 new G4PVPlacement(rotation, 660 G4ThreeVector(0.5 * DWCal1x, -0.5 * DWCal1y, 661 (0.5 * DWCthickness - 0.5 * DWCal1thickness)), 662 DWCal1logical, "DWC_al1_3", fDWClogical, true, 3, true); 663 new G4PVPlacement( 664 0, 665 G4ThreeVector(0.5 * DWCal1y, 0.5 * DWCal1x, 666 -(0.5 * DWCthickness - 0.5 * DWCal1thickness)), 667 DWCal1logical, "DWC_al1_4", fDWClogical, true, 4, true); 668 new G4PVPlacement( 669 rotation, 670 G4ThreeVector(-0.5 * DWCal1x, 0.5 * DWCal1y, 671 -(0.5 * DWCthickness - 0.5 * DWCal1thickness)), 672 DWCal1logical, "DWC_al1_5", fDWClogical, true, 5, true); 673 new G4PVPlacement( 674 0, 675 G4ThreeVector(-0.5 * DWCal1y, -0.5 * DWCal1x, 676 -(0.5 * DWCthickness - 0.5 * DWCal1thickness)), 677 DWCal1logical, "DWC_al1_6", fDWClogical, true, 6, true); 678 new G4PVPlacement( 679 rotation, 680 G4ThreeVector(0.5 * DWCal1x, -0.5 * DWCal1y, 681 -(0.5 * DWCthickness - 0.5 * DWCal1thickness)), 682 DWCal1logical, "DWC_al1_7", fDWClogical, true, 7, true); 683 684 // WChambAl2 685 G4double DWCal2thickness = 2 * 2.25 * CLHEP::cm; 686 G4double DWCal2x = 2 * 10.75 * CLHEP::cm; 687 G4double DWCal2y = 2 * 2.5 * CLHEP::mm; 688 G4Box *DWCal2solid = 689 new G4Box("DWC_al2", 0.5 * DWCal2x, 0.5 * DWCal2y, 0.5 * DWCal2thickness); 690 auto DWCal2logical = new G4LogicalVolume(DWCal2solid, fMatAl, "DWC_al2"); 691 DWCal2logical->SetVisAttributes(visAttributes); 692 new G4PVPlacement(0, G4ThreeVector(0.5 * DWCal2y, 0.5 * DWCal2x, 0), 693 DWCal2logical, "DWC_al2_0", fDWClogical, true, 0, true); 694 new G4PVPlacement(rotation, G4ThreeVector(-0.5 * DWCal2x, 0.5 * DWCal2y, 0), 695 DWCal2logical, "DWC_al2_1", fDWClogical, true, 1, true); 696 new G4PVPlacement(0, G4ThreeVector(-0.5 * DWCal2y, -0.5 * DWCal2x, 0), 697 DWCal2logical, "DWC_al2_2", fDWClogical, true, 2, true); 698 new G4PVPlacement(rotation, G4ThreeVector(0.5 * DWCal2x, -0.5 * DWCal2y, 0), 699 DWCal2logical, "DWC_al2_3", fDWClogical, true, 3, true); 700 701 // WChambGasVet 702 G4double DWCgasVetThickness = 2 * 2.5 * CLHEP::mm; 703 G4double DWCgasVetX = 2 * 7.25 * CLHEP::cm; 704 G4double DWCgasVetY = 2 * 1.25 * CLHEP::cm; 705 G4Box *DWCgasVetSolid = new G4Box("DWC_gasVet", 0.5 * DWCgasVetX, 706 0.5 * DWCgasVetY, 0.5 * DWCgasVetThickness); 707 auto DWCgasVetLogical = 708 new G4LogicalVolume(DWCgasVetSolid, fMatPolyethylene, "DWC_gasVet"); 709 DWCgasVetLogical->SetVisAttributes(visAttributes); 710 new G4PVPlacement(0, 711 G4ThreeVector(0.5 * DWCgasVetY, 0.5 * DWCgasVetX, 712 0.5 * (DWCgasThickness - DWCgasVetThickness)), 713 DWCgasVetLogical, "DWC_gasVet_0", fDWCgasLogical, true, 0, 714 true); 715 new G4PVPlacement(rotation, 716 G4ThreeVector(-0.5 * DWCgasVetX, 0.5 * DWCgasVetY, 717 0.5 * (DWCgasThickness - DWCgasVetThickness)), 718 DWCgasVetLogical, "DWC_gasVet_1", fDWCgasLogical, true, 1, 719 true); 720 new G4PVPlacement(0, 721 G4ThreeVector(-0.5 * DWCgasVetY, -0.5 * DWCgasVetX, 722 0.5 * (DWCgasThickness - DWCgasVetThickness)), 723 DWCgasVetLogical, "DWC_gasVet_2", fDWCgasLogical, true, 2, 724 true); 725 new G4PVPlacement(rotation, 726 G4ThreeVector(0.5 * DWCgasVetX, -0.5 * DWCgasVetY, 727 0.5 * (DWCgasThickness - DWCgasVetThickness)), 728 DWCgasVetLogical, "DWC_gasVet_3", fDWCgasLogical, true, 3, 729 true); 730 new G4PVPlacement(0, G4ThreeVector(0.5 * DWCgasVetY, 0.5 * DWCgasVetX, 0), 731 DWCgasVetLogical, "DWC_gasVet_4", fDWCgasLogical, true, 4, 732 true); 733 new G4PVPlacement( 734 rotation, G4ThreeVector(-0.5 * DWCgasVetX, 0.5 * DWCgasVetY, 0), 735 DWCgasVetLogical, "DWC_gasVet_5", fDWCgasLogical, true, 5, true); 736 new G4PVPlacement(0, G4ThreeVector(-0.5 * DWCgasVetY, -0.5 * DWCgasVetX, 0), 737 DWCgasVetLogical, "DWC_gasVet_6", fDWCgasLogical, true, 6, 738 true); 739 new G4PVPlacement( 740 rotation, G4ThreeVector(0.5 * DWCgasVetX, -0.5 * DWCgasVetY, 0), 741 DWCgasVetLogical, "DWC_gasVet_7", fDWCgasLogical, true, 7, true); 742 new G4PVPlacement( 743 0, 744 G4ThreeVector(0.5 * DWCgasVetY, 0.5 * DWCgasVetX, 745 -0.5 * (DWCgasThickness - DWCgasVetThickness)), 746 DWCgasVetLogical, "DWC_gasVet_8", fDWCgasLogical, true, 8, true); 747 new G4PVPlacement( 748 rotation, 749 G4ThreeVector(-0.5 * DWCgasVetX, 0.5 * DWCgasVetY, 750 -0.5 * (DWCgasThickness - DWCgasVetThickness)), 751 DWCgasVetLogical, "DWC_gasVet_9", fDWCgasLogical, true, 9, true); 752 new G4PVPlacement( 753 0, 754 G4ThreeVector(-0.5 * DWCgasVetY, -0.5 * DWCgasVetX, 755 -0.5 * (DWCgasThickness - DWCgasVetThickness)), 756 DWCgasVetLogical, "DWC_gasVet_10", fDWCgasLogical, true, 10, true); 757 new G4PVPlacement( 758 rotation, 759 G4ThreeVector(0.5 * DWCgasVetX, -0.5 * DWCgasVetY, 760 -0.5 * (DWCgasThickness - DWCgasVetThickness)), 761 DWCgasVetLogical, "DWC_gasVet_11", fDWCgasLogical, true, 11, true); 762 } 763 764 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 765 766 void HGCalTBMaterials::PlaceItemInLogicalVolume(std::string aName, 767 G4double &aZ0, 768 G4LogicalVolume *aLogicMother) { 769 if (aName.find("_DAISY") != std::string::npos) { 770 aName.resize(aName.find("_DAISY")); 771 if (fCopyCounterMap.find(aName) == fCopyCounterMap.end()) 772 fCopyCounterMap[aName] = 0; 773 double dx_ = 2 * std::sin(fAlpha) * 11 * fSiPixelSideLength; 774 double dy_ = 11 * fSiPixelSideLength * (2. + 2 * std::cos(fAlpha)); 775 int nRows_[3] = {1, 2, 1}; 776 for (int nC = 0; nC < 3; nC++) { 777 for (int middle_index = 0; middle_index < nRows_[nC]; middle_index++) { 778 new G4PVPlacement( 779 0, 780 G4ThreeVector(dy_ * (middle_index - nRows_[nC] / 2. + 0.5), 781 -nC * dx_ / 2, aZ0 + 0.5 * fThicknessMap[aName]), 782 fLogicalVolumeMap[aName], aName, aLogicMother, false, 783 fCopyCounterMap[aName]++, true); 784 if (nC <= 0) 785 continue; 786 new G4PVPlacement( 787 0, 788 G4ThreeVector(dy_ * (middle_index - nRows_[nC] / 2. + 0.5), 789 +nC * dx_ / 2, aZ0 + 0.5 * fThicknessMap[aName]), 790 fLogicalVolumeMap[aName], aName, aLogicMother, false, 791 fCopyCounterMap[aName]++, true); 792 } 793 } 794 aZ0 += fThicknessMap[aName]; 795 } else if (aName.find("_SUMMER2017TRIPLET") != std::string::npos) { 796 aName.resize(aName.find("_SUMMER2017TRIPLET")); 797 if (fCopyCounterMap.find(aName) == fCopyCounterMap.end()) 798 fCopyCounterMap[aName] = 0; 799 double dx_ = 2 * std::sin(fAlpha) * 11 * fSiPixelSideLength; 800 double dy_ = 11 * fSiPixelSideLength * (2. + 2 * std::cos(fAlpha)); 801 int nRows_[2] = {1, 2}; 802 for (int nC = 0; nC < 2; nC++) { 803 new G4PVPlacement(0, 804 G4ThreeVector(-dy_ * (0 - nRows_[nC] / 2. + 0.5), 805 -nC * dx_ / 2, 806 aZ0 + 0.5 * fThicknessMap[aName]), 807 fLogicalVolumeMap[aName], aName, aLogicMother, false, 808 fCopyCounterMap[aName]++, true); 809 if (nC <= 0) 810 continue; 811 new G4PVPlacement(0, 812 G4ThreeVector(-dy_ * (0 - nRows_[nC] / 2. + 0.5), 813 +nC * dx_ / 2, 814 aZ0 + 0.5 * fThicknessMap[aName]), 815 fLogicalVolumeMap[aName], aName, aLogicMother, false, 816 fCopyCounterMap[aName]++, true); 817 } 818 aZ0 += fThicknessMap[aName]; 819 } else if (aName.find("_rot30") != std::string::npos) { 820 aName.resize(aName.find("_rot30")); 821 if (fCopyCounterMap.find(aName) == fCopyCounterMap.end()) 822 fCopyCounterMap[aName] = 0; 823 G4RotationMatrix *rot = new G4RotationMatrix; 824 rot->rotateZ(30 * CLHEP::deg); 825 new G4PVPlacement(rot, 826 G4ThreeVector(0, 0, aZ0 + 0.5 * fThicknessMap[aName]), 827 fLogicalVolumeMap[aName], aName, aLogicMother, false, 828 fCopyCounterMap[aName]++, true); 829 aZ0 += fThicknessMap[aName]; 830 } else { 831 if (fCopyCounterMap.find(aName) == fCopyCounterMap.end()) 832 fCopyCounterMap[aName] = 0; 833 new G4PVPlacement(0, G4ThreeVector(0, 0, aZ0 + 0.5 * fThicknessMap[aName]), 834 fLogicalVolumeMap[aName], aName, aLogicMother, false, 835 fCopyCounterMap[aName]++, true); 836 aZ0 += fThicknessMap[aName]; 837 } 838 } 839 840 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 841