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/WLSPrimaryGeneratorAction.cc 28 /// \brief Implementation of the WLSPrimaryGeneratorAction class 29 // 30 // 31 32 #include "WLSPrimaryGeneratorAction.hh" 33 34 #include "WLSDetectorConstruction.hh" 35 #include "WLSPrimaryGeneratorMessenger.hh" 36 37 #include "G4AutoLock.hh" 38 #include "G4Event.hh" 39 #include "G4GeneralParticleSource.hh" 40 #include "G4Material.hh" 41 #include "G4MaterialPropertiesTable.hh" 42 #include "G4OpticalPhoton.hh" 43 #include "G4PhysicsTable.hh" 44 #include "G4SystemOfUnits.hh" 45 #include "G4UImanager.hh" 46 #include "Randomize.hh" 47 48 namespace 49 { 50 G4Mutex gen_mutex = G4MUTEX_INITIALIZER; 51 } 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 G4bool WLSPrimaryGeneratorAction::fFirst = false; 56 57 WLSPrimaryGeneratorAction::WLSPrimaryGeneratorAction(WLSDetectorConstruction* dc) 58 { 59 fDetector = dc; 60 61 fParticleGun = new G4GeneralParticleSource(); 62 fGunMessenger = new WLSPrimaryGeneratorMessenger(this); 63 } 64 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 66 67 WLSPrimaryGeneratorAction::~WLSPrimaryGeneratorAction() 68 { 69 delete fParticleGun; 70 delete fGunMessenger; 71 if (fIntegralTable) { 72 fIntegralTable->clearAndDestroy(); 73 delete fIntegralTable; 74 } 75 } 76 77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 78 79 void WLSPrimaryGeneratorAction::SetDecayTimeConstant(G4double time) 80 { 81 fTimeConstant = time; 82 } 83 84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 85 86 void WLSPrimaryGeneratorAction::BuildEmissionSpectrum() 87 { 88 if (fIntegralTable) return; 89 90 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable(); 91 G4int numOfMaterials = G4Material::GetNumberOfMaterials(); 92 93 if (!fIntegralTable) fIntegralTable = new G4PhysicsTable(numOfMaterials); 94 95 for (G4int i = 0; i < numOfMaterials; ++i) { 96 auto vec = new G4PhysicsFreeVector(); 97 98 G4MaterialPropertiesTable* MPT = (*theMaterialTable)[i]->GetMaterialPropertiesTable(); 99 100 if (MPT) { 101 G4MaterialPropertyVector* theWLSVector = MPT->GetProperty("WLSCOMPONENT"); 102 103 if (theWLSVector) { 104 G4double currentIN = (*theWLSVector)[0]; 105 if (currentIN >= 0.0) { 106 G4double currentPM = theWLSVector->Energy(0); 107 G4double currentCII = 0.0; 108 vec->InsertValues(currentPM, currentCII); 109 G4double prevPM = currentPM; 110 G4double prevCII = currentCII; 111 G4double prevIN = currentIN; 112 113 for (size_t j = 1; j < theWLSVector->GetVectorLength(); ++j) { 114 currentPM = theWLSVector->Energy(j); 115 currentIN = (*theWLSVector)[j]; 116 currentCII = 0.5 * (prevIN + currentIN); 117 currentCII = prevCII + (currentPM - prevPM) * currentCII; 118 vec->InsertValues(currentPM, currentCII); 119 prevPM = currentPM; 120 prevCII = currentCII; 121 prevIN = currentIN; 122 } 123 } 124 } 125 } 126 fIntegralTable->insertAt(i, vec); 127 } 128 } 129 130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 131 132 void WLSPrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent) 133 { 134 if (!fFirst) { 135 fFirst = true; 136 BuildEmissionSpectrum(); 137 } 138 139 if (fUseSampledEnergy) { 140 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable(); 141 142 G4double sampledEnergy = 3. * eV; 143 144 for (size_t j = 0; j < theMaterialTable->size(); ++j) { 145 G4Material* fMaterial = (*theMaterialTable)[j]; 146 if (fMaterial->GetName() == "PMMA") { 147 auto WLSIntensity = fMaterial->GetMaterialPropertiesTable()->GetProperty("WLSCOMPONENT"); 148 149 if (WLSIntensity) { 150 auto WLSIntegral = (G4PhysicsFreeVector*)((*fIntegralTable)(fMaterial->GetIndex())); 151 152 G4double CIImax = WLSIntegral->GetMaxValue(); 153 G4double CIIvalue = G4UniformRand() * CIImax; 154 155 sampledEnergy = WLSIntegral->GetEnergy(CIIvalue); 156 } 157 } 158 } 159 160 // this does not work. 161 G4String cmd = "/gun/energy " + G4UIcommand::ConvertToString(sampledEnergy / eV) + " eV"; 162 G4UImanager::GetUIpointer()->ApplyCommand(cmd); 163 } 164 165 // The code behind this line is not thread safe because polarization 166 // and time are randomly selected and GPS properties are global 167 168 G4AutoLock l(&gen_mutex); 169 if (fParticleGun->GetParticleDefinition() == G4OpticalPhoton::Definition()) { 170 SetOptPhotonPolar(); 171 SetOptPhotonTime(); 172 } 173 174 fParticleGun->GeneratePrimaryVertex(anEvent); 175 } 176 177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 178 179 void WLSPrimaryGeneratorAction::SetOptPhotonPolar() 180 { 181 G4double angle = G4UniformRand() * 360.0 * deg; 182 SetOptPhotonPolar(angle); 183 } 184 185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 186 187 void WLSPrimaryGeneratorAction::SetOptPhotonPolar(G4double angle) 188 { 189 if (fParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") { 190 G4cout << "-> warning from WLSPrimaryGeneratorAction::SetOptPhotonPolar()" 191 << ": the ParticleGun is not an opticalphoton" << G4endl; 192 return; 193 } 194 195 G4ThreeVector normal(1., 0., 0.); 196 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection(); 197 G4ThreeVector product = normal.cross(kphoton); 198 G4double modul2 = product * product; 199 200 G4ThreeVector e_perpend(0., 0., 1.); 201 if (modul2 > 0.) e_perpend = (1. / std::sqrt(modul2)) * product; 202 G4ThreeVector e_paralle = e_perpend.cross(kphoton); 203 204 G4ThreeVector polar = std::cos(angle) * e_paralle + std::sin(angle) * e_perpend; 205 fParticleGun->SetParticlePolarization(polar); 206 } 207 208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 209 210 void WLSPrimaryGeneratorAction::SetOptPhotonTime() 211 { 212 G4double time = -std::log(G4UniformRand()) * fTimeConstant; 213 fParticleGun->SetParticleTime(time); 214 } 215