File indexing completed on 2025-02-23 09:22:25
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032 #include "WLSPrimaryGeneratorAction.hh"
0033
0034 #include "WLSDetectorConstruction.hh"
0035 #include "WLSPrimaryGeneratorMessenger.hh"
0036
0037 #include "G4AutoLock.hh"
0038 #include "G4Event.hh"
0039 #include "G4GeneralParticleSource.hh"
0040 #include "G4Material.hh"
0041 #include "G4MaterialPropertiesTable.hh"
0042 #include "G4OpticalPhoton.hh"
0043 #include "G4PhysicsTable.hh"
0044 #include "G4SystemOfUnits.hh"
0045 #include "G4UImanager.hh"
0046 #include "Randomize.hh"
0047
0048 namespace
0049 {
0050 G4Mutex gen_mutex = G4MUTEX_INITIALIZER;
0051 }
0052
0053
0054
0055 G4bool WLSPrimaryGeneratorAction::fFirst = false;
0056
0057 WLSPrimaryGeneratorAction::WLSPrimaryGeneratorAction(WLSDetectorConstruction* dc)
0058 {
0059 fDetector = dc;
0060
0061 fParticleGun = new G4GeneralParticleSource();
0062 fGunMessenger = new WLSPrimaryGeneratorMessenger(this);
0063 }
0064
0065
0066
0067 WLSPrimaryGeneratorAction::~WLSPrimaryGeneratorAction()
0068 {
0069 delete fParticleGun;
0070 delete fGunMessenger;
0071 if (fIntegralTable) {
0072 fIntegralTable->clearAndDestroy();
0073 delete fIntegralTable;
0074 }
0075 }
0076
0077
0078
0079 void WLSPrimaryGeneratorAction::SetDecayTimeConstant(G4double time)
0080 {
0081 fTimeConstant = time;
0082 }
0083
0084
0085
0086 void WLSPrimaryGeneratorAction::BuildEmissionSpectrum()
0087 {
0088 if (fIntegralTable) return;
0089
0090 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
0091 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
0092
0093 if (!fIntegralTable) fIntegralTable = new G4PhysicsTable(numOfMaterials);
0094
0095 for (G4int i = 0; i < numOfMaterials; ++i) {
0096 auto vec = new G4PhysicsFreeVector();
0097
0098 G4MaterialPropertiesTable* MPT = (*theMaterialTable)[i]->GetMaterialPropertiesTable();
0099
0100 if (MPT) {
0101 G4MaterialPropertyVector* theWLSVector = MPT->GetProperty("WLSCOMPONENT");
0102
0103 if (theWLSVector) {
0104 G4double currentIN = (*theWLSVector)[0];
0105 if (currentIN >= 0.0) {
0106 G4double currentPM = theWLSVector->Energy(0);
0107 G4double currentCII = 0.0;
0108 vec->InsertValues(currentPM, currentCII);
0109 G4double prevPM = currentPM;
0110 G4double prevCII = currentCII;
0111 G4double prevIN = currentIN;
0112
0113 for (size_t j = 1; j < theWLSVector->GetVectorLength(); ++j) {
0114 currentPM = theWLSVector->Energy(j);
0115 currentIN = (*theWLSVector)[j];
0116 currentCII = 0.5 * (prevIN + currentIN);
0117 currentCII = prevCII + (currentPM - prevPM) * currentCII;
0118 vec->InsertValues(currentPM, currentCII);
0119 prevPM = currentPM;
0120 prevCII = currentCII;
0121 prevIN = currentIN;
0122 }
0123 }
0124 }
0125 }
0126 fIntegralTable->insertAt(i, vec);
0127 }
0128 }
0129
0130
0131
0132 void WLSPrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0133 {
0134 if (!fFirst) {
0135 fFirst = true;
0136 BuildEmissionSpectrum();
0137 }
0138
0139 if (fUseSampledEnergy) {
0140 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
0141
0142 G4double sampledEnergy = 3. * eV;
0143
0144 for (size_t j = 0; j < theMaterialTable->size(); ++j) {
0145 G4Material* fMaterial = (*theMaterialTable)[j];
0146 if (fMaterial->GetName() == "PMMA") {
0147 auto WLSIntensity = fMaterial->GetMaterialPropertiesTable()->GetProperty("WLSCOMPONENT");
0148
0149 if (WLSIntensity) {
0150 auto WLSIntegral = (G4PhysicsFreeVector*)((*fIntegralTable)(fMaterial->GetIndex()));
0151
0152 G4double CIImax = WLSIntegral->GetMaxValue();
0153 G4double CIIvalue = G4UniformRand() * CIImax;
0154
0155 sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
0156 }
0157 }
0158 }
0159
0160
0161 G4String cmd = "/gun/energy " + G4UIcommand::ConvertToString(sampledEnergy / eV) + " eV";
0162 G4UImanager::GetUIpointer()->ApplyCommand(cmd);
0163 }
0164
0165
0166
0167
0168 G4AutoLock l(&gen_mutex);
0169 if (fParticleGun->GetParticleDefinition() == G4OpticalPhoton::Definition()) {
0170 SetOptPhotonPolar();
0171 SetOptPhotonTime();
0172 }
0173
0174 fParticleGun->GeneratePrimaryVertex(anEvent);
0175 }
0176
0177
0178
0179 void WLSPrimaryGeneratorAction::SetOptPhotonPolar()
0180 {
0181 G4double angle = G4UniformRand() * 360.0 * deg;
0182 SetOptPhotonPolar(angle);
0183 }
0184
0185
0186
0187 void WLSPrimaryGeneratorAction::SetOptPhotonPolar(G4double angle)
0188 {
0189 if (fParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") {
0190 G4cout << "-> warning from WLSPrimaryGeneratorAction::SetOptPhotonPolar()"
0191 << ": the ParticleGun is not an opticalphoton" << G4endl;
0192 return;
0193 }
0194
0195 G4ThreeVector normal(1., 0., 0.);
0196 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
0197 G4ThreeVector product = normal.cross(kphoton);
0198 G4double modul2 = product * product;
0199
0200 G4ThreeVector e_perpend(0., 0., 1.);
0201 if (modul2 > 0.) e_perpend = (1. / std::sqrt(modul2)) * product;
0202 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
0203
0204 G4ThreeVector polar = std::cos(angle) * e_paralle + std::sin(angle) * e_perpend;
0205 fParticleGun->SetParticlePolarization(polar);
0206 }
0207
0208
0209
0210 void WLSPrimaryGeneratorAction::SetOptPhotonTime()
0211 {
0212 G4double time = -std::log(G4UniformRand()) * fTimeConstant;
0213 fParticleGun->SetParticleTime(time);
0214 }