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 // 27 /// \file optical/wls/src/WLSMaterials.cc 28 /// \brief Implementation of the WLSMaterials class 29 // 30 // 31 #include "WLSMaterials.hh" 32 33 #include "G4NistManager.hh" 34 #include "G4SystemOfUnits.hh" 35 36 WLSMaterials* WLSMaterials::fInstance = nullptr; 37 38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 39 40 WLSMaterials::WLSMaterials() 41 { 42 fNistMan = G4NistManager::Instance(); 43 fNistMan->SetVerbose(2); 44 45 CreateMaterials(); 46 } 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 49 50 WLSMaterials::~WLSMaterials() 51 { 52 delete fAir; 53 delete fPMMA; 54 delete fPethylene; 55 delete fFPethylene; 56 delete fPolystyrene; 57 delete fSilicone; 58 delete fCoating; 59 } 60 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 62 63 WLSMaterials* WLSMaterials::GetInstance() 64 { 65 if (!fInstance) { 66 fInstance = new WLSMaterials(); 67 } 68 return fInstance; 69 } 70 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 72 73 G4Material* WLSMaterials::GetMaterial(const G4String material) 74 { 75 G4Material* mat = fNistMan->FindOrBuildMaterial(material); 76 77 if (!mat) mat = G4Material::GetMaterial(material); 78 if (!mat) { 79 G4ExceptionDescription ed; 80 ed << "Material " << material << " not found!"; 81 G4Exception("WLSMaterials::GetMaterial", "", FatalException, ed); 82 } 83 84 return mat; 85 } 86 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 88 89 void WLSMaterials::CreateMaterials() 90 { 91 G4double density; 92 G4int ncomponents; 93 G4double fractionmass; 94 std::vector<G4int> natoms; 95 std::vector<G4double> fractionMass; 96 std::vector<G4String> elements; 97 98 // Materials Definitions 99 // ===================== 100 101 //-------------------------------------------------- 102 // Vacuum 103 //-------------------------------------------------- 104 105 fNistMan->FindOrBuildMaterial("G4_Galactic"); 106 107 //-------------------------------------------------- 108 // Air 109 //-------------------------------------------------- 110 111 fAir = fNistMan->FindOrBuildMaterial("G4_AIR"); 112 113 //-------------------------------------------------- 114 // WLSfiber PMMA 115 //-------------------------------------------------- 116 117 elements.push_back("C"); 118 natoms.push_back(5); 119 elements.push_back("H"); 120 natoms.push_back(8); 121 elements.push_back("O"); 122 natoms.push_back(2); 123 124 density = 1.190 * g / cm3; 125 126 fPMMA = fNistMan->ConstructNewMaterial("PMMA", elements, natoms, density); 127 128 elements.clear(); 129 natoms.clear(); 130 131 //-------------------------------------------------- 132 // Cladding (polyethylene) 133 //-------------------------------------------------- 134 135 elements.push_back("C"); 136 natoms.push_back(2); 137 elements.push_back("H"); 138 natoms.push_back(4); 139 140 density = 1.200 * g / cm3; 141 142 fPethylene = fNistMan->ConstructNewMaterial("Pethylene", elements, natoms, density); 143 144 elements.clear(); 145 natoms.clear(); 146 147 //-------------------------------------------------- 148 // Double Cladding (fluorinated polyethylene) 149 //-------------------------------------------------- 150 151 elements.push_back("C"); 152 natoms.push_back(2); 153 elements.push_back("H"); 154 natoms.push_back(4); 155 156 density = 1.400 * g / cm3; 157 158 fFPethylene = fNistMan->ConstructNewMaterial("FPethylene", elements, natoms, density); 159 160 elements.clear(); 161 natoms.clear(); 162 163 //-------------------------------------------------- 164 // Polystyrene 165 //-------------------------------------------------- 166 167 elements.push_back("C"); 168 natoms.push_back(8); 169 elements.push_back("H"); 170 natoms.push_back(8); 171 172 density = 1.050 * g / cm3; 173 174 fPolystyrene = fNistMan->ConstructNewMaterial("Polystyrene", elements, natoms, density); 175 176 elements.clear(); 177 natoms.clear(); 178 179 //-------------------------------------------------- 180 // Silicone (Template for Optical Grease) 181 //-------------------------------------------------- 182 183 elements.push_back("C"); 184 natoms.push_back(2); 185 elements.push_back("H"); 186 natoms.push_back(6); 187 188 density = 1.060 * g / cm3; 189 190 fSilicone = fNistMan->ConstructNewMaterial("Silicone", elements, natoms, density); 191 192 elements.clear(); 193 natoms.clear(); 194 195 //-------------------------------------------------- 196 // Aluminium 197 //-------------------------------------------------- 198 199 fNistMan->FindOrBuildMaterial("G4_Al"); 200 201 //-------------------------------------------------- 202 // TiO2 203 //-------------------------------------------------- 204 205 elements.push_back("Ti"); 206 natoms.push_back(1); 207 elements.push_back("O"); 208 natoms.push_back(2); 209 210 density = 4.26 * g / cm3; 211 212 G4Material* TiO2 = fNistMan->ConstructNewMaterial("TiO2", elements, natoms, density); 213 214 elements.clear(); 215 natoms.clear(); 216 217 //-------------------------------------------------- 218 // Scintillator Coating - 15% TiO2 and 85% polystyrene by weight. 219 //-------------------------------------------------- 220 221 density = 1.52 * g / cm3; 222 223 fCoating = new G4Material("Coating", density, ncomponents = 2); 224 225 fCoating->AddMaterial(TiO2, fractionmass = 15 * perCent); 226 fCoating->AddMaterial(fPolystyrene, fractionmass = 85 * perCent); 227 228 // 229 // ------------ Generate & Add Material Properties Table ------------ 230 // 231 232 std::vector<G4double> energy = { 233 2.00 * eV, 2.03 * eV, 2.06 * eV, 2.09 * eV, 2.12 * eV, 2.15 * eV, 2.18 * eV, 2.21 * eV, 234 2.24 * eV, 2.27 * eV, 2.30 * eV, 2.33 * eV, 2.36 * eV, 2.39 * eV, 2.42 * eV, 2.45 * eV, 235 2.48 * eV, 2.51 * eV, 2.54 * eV, 2.57 * eV, 2.60 * eV, 2.63 * eV, 2.66 * eV, 2.69 * eV, 236 2.72 * eV, 2.75 * eV, 2.78 * eV, 2.81 * eV, 2.84 * eV, 2.87 * eV, 2.90 * eV, 2.93 * eV, 237 2.96 * eV, 2.99 * eV, 3.02 * eV, 3.05 * eV, 3.08 * eV, 3.11 * eV, 3.14 * eV, 3.17 * eV, 238 3.20 * eV, 3.23 * eV, 3.26 * eV, 3.29 * eV, 3.32 * eV, 3.35 * eV, 3.38 * eV, 3.41 * eV, 239 3.44 * eV, 3.47 * eV}; 240 241 std::vector<G4double> energySmall = {2.0 * eV, 3.47 * eV}; 242 243 //-------------------------------------------------- 244 // Air 245 //-------------------------------------------------- 246 247 std::vector<G4double> refractiveIndex = {1.0, 1.0}; 248 249 auto mpt = new G4MaterialPropertiesTable(); 250 mpt->AddProperty("RINDEX", energySmall, refractiveIndex); 251 252 fAir->SetMaterialPropertiesTable(mpt); 253 254 //-------------------------------------------------- 255 // PMMA for WLSfibers 256 //-------------------------------------------------- 257 258 std::vector<G4double> refractiveIndexWLSfiber = {1.60, 1.60}; 259 260 std::vector<G4double> absWLSfiber = { 261 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 262 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 263 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 5.40 * m, 264 5.40 * m, 5.40 * m, 1.10 * m, 1.10 * m, 1.10 * m, 1.10 * m, 1.10 * m, 1.10 * m, 1.10 * m, 265 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm, 266 1. * mm, 1. * mm, 1. * mm, 1. * mm, 1. * mm}; 267 268 std::vector<G4double> emissionFib = {0.05, 0.10, 0.30, 0.50, 0.75, 1.00, 1.50, 1.85, 2.30, 2.75, 269 3.25, 3.80, 4.50, 5.20, 6.00, 7.00, 8.50, 9.50, 11.1, 12.4, 270 12.9, 13.0, 12.8, 12.3, 11.1, 11.0, 12.0, 11.0, 17.0, 16.9, 271 15.0, 9.00, 2.50, 1.00, 0.05, 0.00, 0.00, 0.00, 0.00, 0.00, 272 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00}; 273 274 // Add entries into properties table 275 auto mptWLSfiber = new G4MaterialPropertiesTable(); 276 mptWLSfiber->AddProperty("RINDEX", energySmall, refractiveIndexWLSfiber); 277 mptWLSfiber->AddProperty("WLSABSLENGTH", energy, absWLSfiber); 278 mptWLSfiber->AddProperty("WLSCOMPONENT", energy, emissionFib); 279 mptWLSfiber->AddConstProperty("WLSTIMECONSTANT", 0.5 * ns); 280 281 fPMMA->SetMaterialPropertiesTable(mptWLSfiber); 282 283 //-------------------------------------------------- 284 // Polyethylene 285 //-------------------------------------------------- 286 287 std::vector<G4double> refractiveIndexClad1 = {1.49, 1.49}; 288 289 std::vector<G4double> absClad = {20.0 * m, 20.0 * m}; 290 291 // Add entries into properties table 292 auto mptClad1 = new G4MaterialPropertiesTable(); 293 mptClad1->AddProperty("RINDEX", energySmall, refractiveIndexClad1); 294 mptClad1->AddProperty("ABSLENGTH", energySmall, absClad); 295 296 fPethylene->SetMaterialPropertiesTable(mptClad1); 297 298 //-------------------------------------------------- 299 // Fluorinated Polyethylene 300 //-------------------------------------------------- 301 302 std::vector<G4double> refractiveIndexClad2 = {1.42, 1.42}; 303 304 // Add entries into properties table 305 auto mptClad2 = new G4MaterialPropertiesTable(); 306 mptClad2->AddProperty("RINDEX", energySmall, refractiveIndexClad2); 307 mptClad2->AddProperty("ABSLENGTH", energySmall, absClad); 308 309 fFPethylene->SetMaterialPropertiesTable(mptClad2); 310 311 //-------------------------------------------------- 312 // Silicone 313 //-------------------------------------------------- 314 315 std::vector<G4double> refractiveIndexSilicone = {1.46, 1.46}; 316 317 // Add entries into properties table 318 auto mptSilicone = new G4MaterialPropertiesTable(); 319 mptSilicone->AddProperty("RINDEX", energySmall, refractiveIndexSilicone); 320 mptSilicone->AddProperty("ABSLENGTH", energySmall, absClad); 321 322 fSilicone->SetMaterialPropertiesTable(mptSilicone); 323 324 //-------------------------------------------------- 325 // Polystyrene 326 //-------------------------------------------------- 327 328 std::vector<G4double> refractiveIndexPS = {1.50, 1.50}; 329 330 std::vector<G4double> absPS = {2. * cm, 2. * cm}; 331 332 std::vector<G4double> scintilFast = { 333 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 334 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 1.0, 1.0, 1.0, 335 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0}; 336 337 // Add entries into properties table 338 auto mptPolystyrene = new G4MaterialPropertiesTable(); 339 mptPolystyrene->AddProperty("RINDEX", energySmall, refractiveIndexPS); 340 mptPolystyrene->AddProperty("ABSLENGTH", energySmall, absPS); 341 mptPolystyrene->AddProperty("SCINTILLATIONCOMPONENT1", energy, scintilFast); 342 mptPolystyrene->AddConstProperty("SCINTILLATIONYIELD", 10. / keV); 343 mptPolystyrene->AddConstProperty("RESOLUTIONSCALE", 1.0); 344 mptPolystyrene->AddConstProperty("SCINTILLATIONTIMECONSTANT1", 10. * ns); 345 346 fPolystyrene->SetMaterialPropertiesTable(mptPolystyrene); 347 348 // Set the Birks Constant for the Polystyrene scintillator 349 fPolystyrene->GetIonisation()->SetBirksConstant(0.126 * mm / MeV); 350 } 351