Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
WLSPrimaryGeneratorAction Class Reference

#include <WLSPrimaryGeneratorAction.hh>

Inheritance diagram for WLSPrimaryGeneratorAction:
Collaboration diagram for WLSPrimaryGeneratorAction:

Public Member Functions

 WLSPrimaryGeneratorAction (WLSDetectorConstruction *)
 
virtual ~WLSPrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *)
 
void BuildEmissionSpectrum ()
 
void SetOptPhotonPolar (G4double)
 
void SetDecayTimeConstant (G4double)
 
- Public Member Functions inherited from G4VUserPrimaryGeneratorAction
 G4VUserPrimaryGeneratorAction ()
 
virtual ~G4VUserPrimaryGeneratorAction ()
 

Protected Attributes

G4PhysicsTablefIntegralTable
 

Detailed Description

Definition at line 51 of file WLSPrimaryGeneratorAction.hh.

Constructor & Destructor Documentation

WLSPrimaryGeneratorAction::WLSPrimaryGeneratorAction ( WLSDetectorConstruction dc)

Definition at line 65 of file WLSPrimaryGeneratorAction.cc.

66 {
67  fDetector = dc;
68  fIntegralTable = NULL;
69 
70  fParticleGun = new G4GeneralParticleSource();
71 
72  fGunMessenger = new WLSPrimaryGeneratorMessenger(this);
73 
74  // G4String particleName;
75  // G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
76 
77  fTimeConstant = 0.;
78 
79 // fParticleGun->SetParticleDefinition(particleTable->
80 // FindParticle(particleName="opticalphoton"));
81 }
WLSPrimaryGeneratorAction::~WLSPrimaryGeneratorAction ( )
virtual

Definition at line 85 of file WLSPrimaryGeneratorAction.cc.

86 {
87  delete fParticleGun;
88  delete fGunMessenger;
89  if (fIntegralTable) {
91  delete fIntegralTable;
92  }
93 }
void clearAndDestroy()

Here is the call graph for this function:

Member Function Documentation

void WLSPrimaryGeneratorAction::BuildEmissionSpectrum ( )

Definition at line 104 of file WLSPrimaryGeneratorAction.cc.

105 {
106  if (fIntegralTable) return;
107 
108  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
109 
110  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
111 
112  if(!fIntegralTable)fIntegralTable = new G4PhysicsTable(numOfMaterials);
113 
114  for (G4int i=0 ; i < numOfMaterials; i++) {
115 
116  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
118 
119  G4Material* aMaterial = (*theMaterialTable)[i];
120 
121  G4MaterialPropertiesTable* aMaterialPropertiesTable =
122  aMaterial->GetMaterialPropertiesTable();
123 
124  if (aMaterialPropertiesTable) {
125  G4MaterialPropertyVector* theWLSVector =
126  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
127 
128  if (theWLSVector) {
129  G4double currentIN = (*theWLSVector)[0];
130  if (currentIN >= 0.0) {
131  G4double currentPM = theWLSVector->Energy(0);
132  G4double currentCII = 0.0;
133  aPhysicsOrderedFreeVector->
134  InsertValues(currentPM , currentCII);
135  G4double prevPM = currentPM;
136  G4double prevCII = currentCII;
137  G4double prevIN = currentIN;
138 
139  for (size_t j = 1;
140  j < theWLSVector->GetVectorLength();
141  j++)
142  {
143  currentPM = theWLSVector->Energy(j);
144  currentIN = (*theWLSVector)[j];
145  currentCII = 0.5 * (prevIN + currentIN);
146  currentCII = prevCII + (currentPM - prevPM) * currentCII;
147  aPhysicsOrderedFreeVector->
148  InsertValues(currentPM, currentCII);
149  prevPM = currentPM;
150  prevCII = currentCII;
151  prevIN = currentIN;
152  }
153  }
154  }
155  }
156  fIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
157  }
158 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
G4double Energy(size_t index) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void insertAt(size_t, G4PhysicsVector *)
double G4double
Definition: G4Types.hh:76
G4MaterialPropertyVector * GetProperty(const char *key)

Here is the call graph for this function:

Here is the caller graph for this function:

void WLSPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 162 of file WLSPrimaryGeneratorAction.cc.

163 {
164  if (!fFirst) {
165  fFirst = true;
167  }
168 
169 #ifdef use_sampledEnergy
170  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
171 
172  G4double sampledEnergy = 3*eV;
173 
174  for (size_t j=0 ; j<theMaterialTable->size() ; j++) {
175  G4Material* fMaterial = (*theMaterialTable)[j];
176  if (fMaterial->GetName() == "PMMA" ) {
177  G4MaterialPropertiesTable* aMaterialPropertiesTable =
178  fMaterial->GetMaterialPropertiesTable();
179  const G4MaterialPropertyVector* WLSIntensity =
180  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
181 
182  if (WLSIntensity) {
183  G4int MaterialIndex = fMaterial->GetIndex();
184  G4PhysicsOrderedFreeVector* WLSIntegral =
185  (G4PhysicsOrderedFreeVector*)((*fIntegralTable)(MaterialIndex));
186 
187  G4double CIImax = WLSIntegral->GetMaxValue();
188  G4double CIIvalue = G4UniformRand()*CIImax;
189 
190  sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
191  }
192  }
193  }
194 
195  //fParticleGun->SetParticleEnergy(sampledEnergy);
196 #endif
197 
198  //The code behing this line is not thread safe because polarization
199  //and time are randomly selected and GPS properties are global
200  G4AutoLock l(&gen_mutex);
201  if(fParticleGun->GetParticleDefinition()->GetParticleName()=="opticalphoton"){
203  SetOptPhotonTime();
204  }
205 
206  fParticleGun->GeneratePrimaryVertex(anEvent);
207 }
size_t GetIndex() const
Definition: G4Material.hh:262
const G4String & GetName() const
Definition: G4Material.hh:178
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
#define G4UniformRand()
Definition: Randomize.hh:97
static constexpr double eV
Definition: G4SIunits.hh:215
G4ParticleDefinition * GetParticleDefinition() const
G4double GetEnergy(G4double aValue)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
double G4double
Definition: G4Types.hh:76
G4MaterialPropertyVector * GetProperty(const char *key)

Here is the call graph for this function:

void WLSPrimaryGeneratorAction::SetDecayTimeConstant ( G4double  time)

Definition at line 97 of file WLSPrimaryGeneratorAction.cc.

98 {
99  fTimeConstant = time;
100 }
void WLSPrimaryGeneratorAction::SetOptPhotonPolar ( G4double  angle)

Definition at line 219 of file WLSPrimaryGeneratorAction.cc.

220 {
221  if (fParticleGun->GetParticleDefinition()->GetParticleName()!="opticalphoton")
222  {
223  G4cout << "-> warning from WLSPrimaryGeneratorAction::SetOptPhotonPolar()"
224  << ": the ParticleGun is not an opticalphoton" << G4endl;
225  return;
226  }
227 
228  G4ThreeVector normal (1., 0., 0.);
229  G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
230  G4ThreeVector product = normal.cross(kphoton);
231  G4double modul2 = product*product;
232 
233  G4ThreeVector e_perpend (0., 0., 1.);
234  if (modul2 > 0.) e_perpend = (1./std::sqrt(modul2))*product;
235  G4ThreeVector e_paralle = e_perpend.cross(kphoton);
236 
237  G4ThreeVector polar = std::cos(angle)*e_paralle + std::sin(angle)*e_perpend;
238  fParticleGun->SetParticlePolarization(polar);
239 
240 }
static G4double angle[DIM]
const G4String & GetParticleName() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4GLOB_DLL std::ostream G4cout
void SetParticlePolarization(G4ThreeVector aVal)
G4ParticleDefinition * GetParticleDefinition() const
G4ThreeVector GetParticleMomentumDirection() const
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4PhysicsTable* WLSPrimaryGeneratorAction::fIntegralTable
protected

Definition at line 70 of file WLSPrimaryGeneratorAction.hh.


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