Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: WLSPrimaryGeneratorAction.cc 69561 2013-05-08 12:25:56Z gcosmo $
27 //
28 /// \file optical/wls/src/WLSPrimaryGeneratorAction.cc
29 /// \brief Implementation of the WLSPrimaryGeneratorAction class
30 //
31 //
32 #include "G4ios.hh"
33 #include "G4Event.hh"
34 
36 
37 #include "G4Material.hh"
39 
40 #include "G4ParticleTable.hh"
41 #include "G4ParticleDefinition.hh"
42 
43 #include "G4PhysicsTable.hh"
44 
45 #include "Randomize.hh"
46 
48 
51 
52 #include "G4SystemOfUnits.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 G4bool WLSPrimaryGeneratorAction::fFirst = false;
57 
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 }
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78 
80 {
81  delete fParticleGun;
82  delete fGunMessenger;
83  if (fIntegralTable) {
85  delete fIntegralTable;
86  }
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90 
92 {
93  fTimeConstant = time;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97 
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 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
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 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202 
203 void WLSPrimaryGeneratorAction::SetOptPhotonPolar()
204 {
205  G4double angle = G4UniformRand() * 360.0*deg;
206  SetOptPhotonPolar(angle);
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
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 }
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235 
236 void WLSPrimaryGeneratorAction::SetOptPhotonTime()
237 {
238  G4double time = -std::log(G4UniformRand())*fTimeConstant;
239  fParticleGun->SetParticleTime(time);
240 }
G4ParticleDefinition * GetParticleDefinition()
G4MaterialPropertyVector * GetProperty(const char *key)
Definition of the WLSPrimaryGeneratorAction class.
size_t GetIndex() const
Definition: G4Material.hh:260
const G4String & GetName() const
Definition: G4Material.hh:176
Definition of the WLSPrimaryGeneratorMessenger class.
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
Definition of the WLSDetectorConstruction class.
void SetParticleTime(G4double aTime)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetParticlePolarization(G4ThreeVector aVal)
G4double Energy(size_t index) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:571
G4double GetEnergy(G4double aValue)
WLSPrimaryGeneratorAction(WLSDetectorConstruction *)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:250
G4ThreeVector GetParticleMomentumDirection()
virtual void GeneratePrimaries(G4Event *)
void insertAt(size_t, G4PhysicsVector *)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void clearAndDestroy()