Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
WLSPrimaryGeneratorAction Class Reference

#include <Doxymodules_optical.h>

Inheritance diagram for WLSPrimaryGeneratorAction:
G4VUserPrimaryGeneratorAction

Public Member Functions

 WLSPrimaryGeneratorAction (WLSDetectorConstruction *)
 
 ~WLSPrimaryGeneratorAction () override
 
void GeneratePrimaries (G4Event *) override
 
void BuildEmissionSpectrum ()
 
void SetOptPhotonPolar (G4double)
 
void SetDecayTimeConstant (G4double)
 
void SetUseSampledEnergy (G4bool v)
 

Protected Attributes

G4PhysicsTablefIntegralTable = nullptr
 

Private Member Functions

void SetOptPhotonPolar ()
 
void SetOptPhotonTime ()
 

Private Attributes

WLSDetectorConstructionfDetector = nullptr
 
G4GeneralParticleSourcefParticleGun = nullptr
 
WLSPrimaryGeneratorMessengerfGunMessenger = nullptr
 
G4double fTimeConstant = 0.
 
G4bool fUseSampledEnergy = false
 

Static Private Attributes

static G4bool fFirst = false
 

Detailed Description

Definition at line 96 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ WLSPrimaryGeneratorAction()

WLSPrimaryGeneratorAction::WLSPrimaryGeneratorAction ( WLSDetectorConstruction dc)

◆ ~WLSPrimaryGeneratorAction()

WLSPrimaryGeneratorAction::~WLSPrimaryGeneratorAction ( )
override

Definition at line 68 of file WLSPrimaryGeneratorAction.cc.

69{
70 delete fParticleGun;
71 delete fGunMessenger;
73 {
74 fIntegralTable->clearAndDestroy();
75 delete fIntegralTable;
76 }
77}

Member Function Documentation

◆ GeneratePrimaries()

void WLSPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
override

Definition at line 142 of file WLSPrimaryGeneratorAction.cc.

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}

◆ BuildEmissionSpectrum()

void WLSPrimaryGeneratorAction::BuildEmissionSpectrum ( )

Definition at line 88 of file WLSPrimaryGeneratorAction.cc.

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}

◆ SetOptPhotonPolar() [1/2]

void WLSPrimaryGeneratorAction::SetOptPhotonPolar ( G4double  angle)

Definition at line 205 of file WLSPrimaryGeneratorAction.cc.

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}

◆ SetDecayTimeConstant()

void WLSPrimaryGeneratorAction::SetDecayTimeConstant ( G4double  time)

Definition at line 81 of file WLSPrimaryGeneratorAction.cc.

◆ SetUseSampledEnergy()

void WLSPrimaryGeneratorAction::SetUseSampledEnergy ( G4bool  v)
inline

Definition at line 59 of file WLSPrimaryGeneratorAction.hh.

59{ fUseSampledEnergy = v; }

◆ SetOptPhotonPolar() [2/2]

void WLSPrimaryGeneratorAction::SetOptPhotonPolar ( )
private

Definition at line 197 of file WLSPrimaryGeneratorAction.cc.

198{
199 G4double angle = G4UniformRand() * 360.0 * deg;
200 SetOptPhotonPolar(angle);
201}

◆ SetOptPhotonTime()

void WLSPrimaryGeneratorAction::SetOptPhotonTime ( )
private

Definition at line 232 of file WLSPrimaryGeneratorAction.cc.

233{
234 G4double time = -std::log(G4UniformRand()) * fTimeConstant;
235 fParticleGun->SetParticleTime(time);
236}

Member Data Documentation

◆ fIntegralTable

G4PhysicsTable* WLSPrimaryGeneratorAction::fIntegralTable = nullptr
protected

Definition at line 62 of file WLSPrimaryGeneratorAction.hh.

◆ fDetector

WLSDetectorConstruction* WLSPrimaryGeneratorAction::fDetector = nullptr
private

Definition at line 68 of file WLSPrimaryGeneratorAction.hh.

◆ fParticleGun

G4GeneralParticleSource* WLSPrimaryGeneratorAction::fParticleGun = nullptr
private

Definition at line 69 of file WLSPrimaryGeneratorAction.hh.

◆ fGunMessenger

WLSPrimaryGeneratorMessenger* WLSPrimaryGeneratorAction::fGunMessenger = nullptr
private

Definition at line 70 of file WLSPrimaryGeneratorAction.hh.

◆ fFirst

G4bool WLSPrimaryGeneratorAction::fFirst = false
staticprivate

Definition at line 72 of file WLSPrimaryGeneratorAction.hh.

◆ fTimeConstant

G4double WLSPrimaryGeneratorAction::fTimeConstant = 0.
private

Definition at line 73 of file WLSPrimaryGeneratorAction.hh.

◆ fUseSampledEnergy

G4bool WLSPrimaryGeneratorAction::fUseSampledEnergy = false
private

Definition at line 74 of file WLSPrimaryGeneratorAction.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration