Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes
WLSPrimaryGeneratorAction Class Reference

#include <WLSPrimaryGeneratorAction.hh>

Inheritance diagram for WLSPrimaryGeneratorAction:
G4VUserPrimaryGeneratorAction

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 59 of file WLSPrimaryGeneratorAction.cc.

References fIntegralTable.

60 {
61  fDetector = dc;
62  fIntegralTable = NULL;
63 
64  fParticleGun = new G4GeneralParticleSource();
65 
66  fGunMessenger = new WLSPrimaryGeneratorMessenger(this);
67 
68  // G4String particleName;
69  // G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
70 
71  fTimeConstant = 0.;
72 
73 // fParticleGun->SetParticleDefinition(particleTable->
74 // FindParticle(particleName="opticalphoton"));
75 }
WLSPrimaryGeneratorAction::~WLSPrimaryGeneratorAction ( )
virtual

Definition at line 79 of file WLSPrimaryGeneratorAction.cc.

References G4PhysicsTable::clearAndDestroy(), and fIntegralTable.

80 {
81  delete fParticleGun;
82  delete fGunMessenger;
83  if (fIntegralTable) {
85  delete fIntegralTable;
86  }
87 }
void clearAndDestroy()

Member Function Documentation

void WLSPrimaryGeneratorAction::BuildEmissionSpectrum ( )

Definition at line 98 of file WLSPrimaryGeneratorAction.cc.

References G4PhysicsVector::Energy(), fIntegralTable, G4Material::GetMaterialPropertiesTable(), G4Material::GetMaterialTable(), G4Material::GetNumberOfMaterials(), G4MaterialPropertiesTable::GetProperty(), G4PhysicsVector::GetVectorLength(), and G4PhysicsTable::insertAt().

Referenced by GeneratePrimaries().

99 {
100  if (fIntegralTable) return;
101 
102  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
103 
104  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
105 
106  if(!fIntegralTable)fIntegralTable = new G4PhysicsTable(numOfMaterials);
107 
108  for (G4int i=0 ; i < numOfMaterials; i++) {
109 
110  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
112 
113  G4Material* aMaterial = (*theMaterialTable)[i];
114 
115  G4MaterialPropertiesTable* aMaterialPropertiesTable =
116  aMaterial->GetMaterialPropertiesTable();
117 
118  if (aMaterialPropertiesTable) {
119  G4MaterialPropertyVector* theWLSVector =
120  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
121 
122  if (theWLSVector) {
123  G4double currentIN = (*theWLSVector)[0];
124  if (currentIN >= 0.0) {
125  G4double currentPM = theWLSVector->Energy(0);
126  G4double currentCII = 0.0;
127  aPhysicsOrderedFreeVector->
128  InsertValues(currentPM , currentCII);
129  G4double prevPM = currentPM;
130  G4double prevCII = currentCII;
131  G4double prevIN = currentIN;
132 
133  for (size_t j = 1;
134  j < theWLSVector->GetVectorLength();
135  j++)
136  {
137  currentPM = theWLSVector->Energy(j);
138  currentIN = (*theWLSVector)[j];
139  currentCII = 0.5 * (prevIN + currentIN);
140  currentCII = prevCII + (currentPM - prevPM) * currentCII;
141  aPhysicsOrderedFreeVector->
142  InsertValues(currentPM, currentCII);
143  prevPM = currentPM;
144  prevCII = currentCII;
145  prevIN = currentIN;
146  }
147  }
148  }
149  }
150  fIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
151  }
152 }
G4MaterialPropertyVector * GetProperty(const char *key)
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
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:571
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:250
void insertAt(size_t, G4PhysicsVector *)
double G4double
Definition: G4Types.hh:76
void WLSPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 156 of file WLSPrimaryGeneratorAction.cc.

References BuildEmissionSpectrum(), python.hepunit::eV, fIntegralTable, G4UniformRand, G4GeneralParticleSource::GeneratePrimaryVertex(), G4PhysicsOrderedFreeVector::GetEnergy(), G4Material::GetIndex(), G4Material::GetMaterialPropertiesTable(), G4Material::GetMaterialTable(), G4PhysicsOrderedFreeVector::GetMaxValue(), G4Material::GetName(), G4GeneralParticleSource::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4MaterialPropertiesTable::GetProperty(), and SetOptPhotonPolar().

157 {
158  if (!fFirst) {
159  fFirst = true;
161  }
162 
163 
164 #ifdef use_sampledEnergy
165  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
166 
167  G4double sampledEnergy = 3*eV;
168 
169  for (size_t j=0 ; j<theMaterialTable->size() ; j++) {
170  G4Material* fMaterial = (*theMaterialTable)[j];
171  if (fMaterial->GetName() == "PMMA" ) {
172  G4MaterialPropertiesTable* aMaterialPropertiesTable =
173  fMaterial->GetMaterialPropertiesTable();
174  const G4MaterialPropertyVector* WLSIntensity =
175  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
176 
177  if (WLSIntensity) {
178  G4int MaterialIndex = fMaterial->GetIndex();
179  G4PhysicsOrderedFreeVector* WLSIntegral =
180  (G4PhysicsOrderedFreeVector*)((*fIntegralTable)(MaterialIndex));
181 
182  G4double CIImax = WLSIntegral->GetMaxValue();
183  G4double CIIvalue = G4UniformRand()*CIImax;
184 
185  sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
186  }
187  }
188  }
189 
190  //fParticleGun->SetParticleEnergy(sampledEnergy);
191 #endif
192 
193  if(fParticleGun->GetParticleDefinition()->GetParticleName()=="opticalphoton"){
195  SetOptPhotonTime();
196  }
197 
198  fParticleGun->GeneratePrimaryVertex(anEvent);
199 }
G4ParticleDefinition * GetParticleDefinition()
G4MaterialPropertyVector * GetProperty(const char *key)
size_t GetIndex() const
Definition: G4Material.hh:260
const G4String & GetName() const
Definition: G4Material.hh:176
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4double GetEnergy(G4double aValue)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:250
double G4double
Definition: G4Types.hh:76
void WLSPrimaryGeneratorAction::SetDecayTimeConstant ( G4double  time)

Definition at line 91 of file WLSPrimaryGeneratorAction.cc.

92 {
93  fTimeConstant = time;
94 }
void WLSPrimaryGeneratorAction::SetOptPhotonPolar ( G4double  angle)

Definition at line 211 of file WLSPrimaryGeneratorAction.cc.

References CLHEP::Hep3Vector::cross(), G4cout, G4endl, G4GeneralParticleSource::GetParticleDefinition(), G4GeneralParticleSource::GetParticleMomentumDirection(), G4ParticleDefinition::GetParticleName(), SetOptPhotonPolar(), and G4GeneralParticleSource::SetParticlePolarization().

Referenced by GeneratePrimaries(), and SetOptPhotonPolar().

212 {
213  if (fParticleGun->GetParticleDefinition()->GetParticleName()!="opticalphoton")
214  {
215  G4cout << "-> warning from WLSPrimaryGeneratorAction::SetOptPhotonPolar()"
216  << ": the ParticleGun is not an opticalphoton" << G4endl;
217  return;
218  }
219 
220  G4ThreeVector normal (1., 0., 0.);
221  G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
222  G4ThreeVector product = normal.cross(kphoton);
223  G4double modul2 = product*product;
224 
225  G4ThreeVector e_perpend (0., 0., 1.);
226  if (modul2 > 0.) e_perpend = (1./std::sqrt(modul2))*product;
227  G4ThreeVector e_paralle = e_perpend.cross(kphoton);
228 
229  G4ThreeVector polar = std::cos(angle)*e_paralle + std::sin(angle)*e_perpend;
230  fParticleGun->SetParticlePolarization(polar);
231 
232 }
G4ParticleDefinition * GetParticleDefinition()
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
void SetParticlePolarization(G4ThreeVector aVal)
G4ThreeVector GetParticleMomentumDirection()
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76

Field Documentation

G4PhysicsTable* WLSPrimaryGeneratorAction::fIntegralTable
protected

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