Loading...
Searching...
No Matches
WLSPrimaryGeneratorAction.cc
Go to the documentation of this file.
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
33
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
48namespace
49{
50 G4Mutex gen_mutex = G4MUTEX_INITIALIZER;
51}
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
56
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
69{
70 delete fParticleGun;
71 delete fGunMessenger;
73 {
74 fIntegralTable->clearAndDestroy();
75 delete fIntegralTable;
76 }
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
89{
91 return;
92
93 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
94 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
95
97 fIntegralTable = new G4PhysicsTable(numOfMaterials);
98
99 for(G4int i = 0; i < numOfMaterials; ++i)
100 {
101 auto vec = new G4PhysicsFreeVector();
102
104 (*theMaterialTable)[i]->GetMaterialPropertiesTable();
105
106 if(MPT)
107 {
108 G4MaterialPropertyVector* theWLSVector = MPT->GetProperty("WLSCOMPONENT");
109
110 if(theWLSVector)
111 {
112 G4double currentIN = (*theWLSVector)[0];
113 if(currentIN >= 0.0)
114 {
115 G4double currentPM = theWLSVector->Energy(0);
116 G4double currentCII = 0.0;
117 vec->InsertValues(currentPM, currentCII);
118 G4double prevPM = currentPM;
119 G4double prevCII = currentCII;
120 G4double prevIN = currentIN;
121
122 for(size_t j = 1; j < theWLSVector->GetVectorLength(); ++j)
123 {
124 currentPM = theWLSVector->Energy(j);
125 currentIN = (*theWLSVector)[j];
126 currentCII = 0.5 * (prevIN + currentIN);
127 currentCII = prevCII + (currentPM - prevPM) * currentCII;
128 vec->InsertValues(currentPM, currentCII);
129 prevPM = currentPM;
130 prevCII = currentCII;
131 prevIN = currentIN;
132 }
133 }
134 }
135 }
136 fIntegralTable->insertAt(i, vec);
137 }
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
143{
144 if(!fFirst)
145 {
146 fFirst = true;
148 }
149
151 {
152 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
153
154 G4double sampledEnergy = 3. * eV;
155
156 for(size_t j = 0; j < theMaterialTable->size(); ++j)
157 {
158 G4Material* fMaterial = (*theMaterialTable)[j];
159 if(fMaterial->GetName() == "PMMA")
160 {
161 auto WLSIntensity =
162 fMaterial->GetMaterialPropertiesTable()->GetProperty("WLSCOMPONENT");
163
164 if(WLSIntensity)
165 {
166 auto WLSIntegral = (G4PhysicsFreeVector*) ((*fIntegralTable)(
167 fMaterial->GetIndex()));
168
169 G4double CIImax = WLSIntegral->GetMaxValue();
170 G4double CIIvalue = G4UniformRand() * CIImax;
171
172 sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
173 }
174 }
175 }
176
177 // this does not work.
178 G4String cmd = "/gun/energy " + G4UIcommand::ConvertToString(sampledEnergy / eV) + " eV";
179 G4UImanager::GetUIpointer()->ApplyCommand(cmd);
180 }
181
182 // The code behind this line is not thread safe because polarization
183 // and time are randomly selected and GPS properties are global
184
185 G4AutoLock l(&gen_mutex);
186 if(fParticleGun->GetParticleDefinition() == G4OpticalPhoton::Definition())
187 {
190 }
191
192 fParticleGun->GeneratePrimaryVertex(anEvent);
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196
198{
199 G4double angle = G4UniformRand() * 360.0 * deg;
200 SetOptPhotonPolar(angle);
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204
206{
207 if(fParticleGun->GetParticleDefinition()->GetParticleName() !=
208 "opticalphoton")
209 {
210 G4cout << "-> warning from WLSPrimaryGeneratorAction::SetOptPhotonPolar()"
211 << ": the ParticleGun is not an opticalphoton" << G4endl;
212 return;
213 }
214
215 G4ThreeVector normal(1., 0., 0.);
216 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
217 G4ThreeVector product = normal.cross(kphoton);
218 G4double modul2 = product * product;
219
220 G4ThreeVector e_perpend(0., 0., 1.);
221 if(modul2 > 0.)
222 e_perpend = (1. / std::sqrt(modul2)) * product;
223 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
224
225 G4ThreeVector polar =
226 std::cos(angle) * e_paralle + std::sin(angle) * e_perpend;
227 fParticleGun->SetParticlePolarization(polar);
228}
229
230//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231
233{
234 G4double time = -std::log(G4UniformRand()) * fTimeConstant;
235 fParticleGun->SetParticleTime(time);
236}
Definition of the WLSDetectorConstruction class.
Definition of the WLSPrimaryGeneratorAction class.
Definition of the WLSPrimaryGeneratorMessenger class.
WLSPrimaryGeneratorMessenger * fGunMessenger
WLSDetectorConstruction * fDetector
void GeneratePrimaries(G4Event *) override
G4GeneralParticleSource * fParticleGun
WLSPrimaryGeneratorAction(WLSDetectorConstruction *)

Applications | User Support | Publications | Collaboration