source: trunk/examples/advanced/xray_fluorescence/src/XrayFluoPlanePrimaryGeneratorAction.cc @ 807

Last change on this file since 807 was 807, checked in by garnier, 16 years ago

update

File size: 7.2 KB
Line 
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// $Id: XrayFluoPlanePrimaryGeneratorAction.cc
28// GEANT4 tag $Name:
29//
30// Author: Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
31//
32// History:
33// -----------
34// 02 Sep 2003 Alfonso Mantero created
35//
36// -------------------------------------------------------------------
37
38#include "XrayFluoPlanePrimaryGeneratorAction.hh"
39#include "G4DataVector.hh"
40#include "XrayFluoPlaneDetectorConstruction.hh"
41#include "XrayFluoPlanePrimaryGeneratorMessenger.hh"
42#include "XrayFluoRunAction.hh"
43#include "G4Event.hh"
44#include "G4ParticleGun.hh"
45#include "G4ParticleTable.hh"
46#include "G4ParticleDefinition.hh"
47#include "Randomize.hh"
48#include "XrayFluoAnalysisManager.hh"
49#include "XrayFluoDataSet.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52
53XrayFluoPlanePrimaryGeneratorAction::XrayFluoPlanePrimaryGeneratorAction(XrayFluoPlaneDetectorConstruction* XrayFluoDC)
54  :rndmFlag("on"),beam("off"),spectrum("off"),isoVert("off")
55{
56
57  XrayFluoDetector = XrayFluoDC;
58
59  G4int n_particle = 1;
60  particleGun  = new G4ParticleGun(n_particle);
61 
62  //create a messenger for this class
63  gunMessenger = new XrayFluoPlanePrimaryGeneratorMessenger(this);
64  runManager = new XrayFluoRunAction();
65 
66  // default particle kinematic
67 
68  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
69  G4String particleName;
70  G4ParticleDefinition* particle
71    = particleTable->FindParticle(particleName="gamma");
72  particleGun->SetParticleDefinition(particle);
73  particleGun->SetParticleMomentumDirection(G4ThreeVector(0.,0.,1.));
74 
75
76  particleGun->SetParticleEnergy(10.*keV);
77  G4double position = -0.5*(XrayFluoDetector->GetWorldSizeZ());
78  particleGun->SetParticlePosition(G4ThreeVector(0.*cm,0.*cm,position));
79
80  G4cout << "XrayFluoPlanePrimaryGeneratorAction created  UUUUUUUUUUAAAAAAAAAAAAAAAAAAAAAAAaa" << G4endl;
81 
82}
83
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
87XrayFluoPlanePrimaryGeneratorAction::~XrayFluoPlanePrimaryGeneratorAction()
88{
89  delete particleGun;
90  delete gunMessenger;
91  delete runManager;
92
93  G4cout << "XrayFluoPlanePrimaryGeneratorAction deleted" << G4endl;
94
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
99void XrayFluoPlanePrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
100{
101  //this function is called at the begining of event
102  //
103  G4double z0 = -0.5*(XrayFluoDetector->GetWorldSizeZ());
104  G4double y0 = 0.*m, x0 = 0.*m;
105  G4double dX = 0.5*(XrayFluoDetector->GetWorldSizeXY())-(XrayFluoDetector->GetPlaneSizeXY());
106  if (rndmFlag == "on")
107
108    {y0 = (XrayFluoDetector->GetPlaneSizeXY())*(G4UniformRand()-0.5); 
109    x0 = (XrayFluoDetector->GetPlaneSizeXY())*(G4UniformRand()-0.5) + dX; 
110    } 
111
112  z0 = -1 * dX;
113
114  particleGun->SetParticleMomentumDirection(G4ThreeVector(-1.,0.,1.));
115
116  particleGun->SetParticlePosition(G4ThreeVector(x0,y0,z0));
117 
118  //randomize starting point
119  if (beam == "on")
120    {
121      G4double radius = 0.5 * mm;
122      G4double rho = radius*std::sqrt(G4UniformRand());
123      G4double theta = 2*pi*G4UniformRand()*rad;
124      G4double position = -0.5*(XrayFluoDetector->GetWorldSizeZ());
125     
126      G4double y = rho * std::sin(theta);
127      G4double x = rho * std::cos(theta);
128     
129      particleGun->SetParticlePosition(G4ThreeVector(x,y,position));
130    }
131  //shoot particles according to a certain spectrum
132  if (spectrum =="on")
133    {
134      G4String particle =  particleGun->GetParticleDefinition()
135        ->GetParticleName();
136      if(particle == "proton"|| particle == "alpha")
137        {
138          G4DataVector* energies =  runManager->GetEnergies();
139          G4DataVector* data =  runManager->GetData();
140         
141          G4double sum = runManager->GetDataSum();
142          G4double partSum = 0;
143          G4int j = 0;
144          G4double random= sum*G4UniformRand();
145          while (partSum<random)
146            {
147              partSum += (*data)[j];
148              j++;
149            }
150         
151          particleGun->SetParticleEnergy((*energies)[j]);
152       
153        }
154      else if (particle == "gamma")
155        {
156          const XrayFluoDataSet* dataSet = runManager->GetGammaSet();
157         
158          G4int i = 0;
159          G4int id = 0;
160          G4double minEnergy = 0. * keV;
161          G4double particleEnergy= 0.;
162          G4double maxEnergy = 10. * keV;
163          G4double energyRange = maxEnergy - minEnergy;
164
165           while ( i == 0)
166            {
167              G4double random = G4UniformRand();
168             
169              G4double randomNum = G4UniformRand(); //*5.0E6;
170             
171              particleEnergy = (random*energyRange) + minEnergy;
172             
173              if ((dataSet->FindValue(particleEnergy,id)) > randomNum)
174                {
175                  i = 1;
176                 
177                }
178            }
179           particleGun->SetParticleEnergy(particleEnergy);
180        }
181    }
182 
183  if (isoVert == "on")
184    {
185      G4double rho = 1. *m;
186      //theta in [0;pi/2]
187      G4double theta = (pi/2)*G4UniformRand();
188      //phi in [-pi;pi]
189      G4double phi = (G4UniformRand()*2*pi)- pi;
190      G4double x = rho*std::sin(theta)*std::sin(phi);
191      G4double y = rho*std::sin(theta)*std::cos(phi);
192      G4double z = -(rho*std::cos(theta));
193      particleGun->SetParticlePosition(G4ThreeVector(x,y,z));
194     
195      G4double Xdim = XrayFluoDetector->GetPlaneSizeXY();
196      G4double Ydim = XrayFluoDetector->GetPlaneSizeXY();
197     
198      G4double Dx = Xdim*(G4UniformRand()-0.5);
199     
200      G4double Dy = Ydim*(G4UniformRand()-0.5);
201     
202      particleGun->SetParticleMomentumDirection(G4ThreeVector(-x+Dx,-y+Dy,-z));
203     
204    }
205#ifdef G4ANALYSIS_USE
206
207  G4double partEnergy = particleGun->GetParticleEnergy();
208  XrayFluoAnalysisManager* analysis =  XrayFluoAnalysisManager::getInstance();
209  analysis->analysePrimaryGenerator(partEnergy/keV);
210
211#endif
212 
213  particleGun->GeneratePrimaryVertex(anEvent);
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217
218
219
220
221
222
223
224
225
Note: See TracBrowser for help on using the repository browser.