source: trunk/source/event/src/G4SingleParticleSource.cc @ 1271

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

import all except CVS

File size: 6.0 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///////////////////////////////////////////////////////////////////////////////
28//
29// MODULE:        G4SingleParticleSource.hh
30//
31// Version:      1.0
32// Date:         5/02/04
33// Author:       Fan Lei
34// Organisation: QinetiQ ltd.
35// Customer:     ESA/ESTEC
36//
37///////////////////////////////////////////////////////////////////////////////
38//
39// CHANGE HISTORY
40// --------------
41//
42// Version 1.0, 05/02/2004, Fan Lei, Created.
43//      Based on the G4GeneralParticleSource class in Geant4 v6.0
44//
45///////////////////////////////////////////////////////////////////////////////
46//
47#include "G4PrimaryParticle.hh"
48#include "G4Event.hh"
49#include "Randomize.hh"
50#include <cmath>
51#include "G4ParticleTable.hh"
52#include "G4Geantino.hh"
53#include "G4ParticleDefinition.hh"
54#include "G4IonTable.hh"
55#include "G4Ions.hh"
56#include "G4TrackingManager.hh"
57#include "G4Track.hh"
58#include "G4SingleParticleSource.hh"
59
60G4SingleParticleSource::G4SingleParticleSource()
61{
62  // Initialise all variables
63  // Position distribution Variables
64
65  NumberOfParticlesToBeGenerated = 1;
66  particle_definition = G4Geantino::GeantinoDefinition();
67  G4ThreeVector zero;
68  particle_momentum_direction = G4ParticleMomentum(1,0,0);
69  particle_energy = 1.0*MeV;
70  particle_position = zero;
71  particle_time = 0.0;
72  particle_polarization = zero;
73  particle_charge = 0.0;
74  particle_weight = 1.0;
75
76  biasRndm = new G4SPSRandomGenerator();
77  posGenerator = new G4SPSPosDistribution();
78  posGenerator->SetBiasRndm(biasRndm);
79  angGenerator = new G4SPSAngDistribution();
80  angGenerator->SetPosDistribution(posGenerator);
81  angGenerator->SetBiasRndm(biasRndm);
82  eneGenerator = new G4SPSEneDistribution();
83  eneGenerator->SetBiasRndm(biasRndm);
84
85  // verbosity
86  verbosityLevel = 0;
87
88}
89
90G4SingleParticleSource::~G4SingleParticleSource()
91{}
92
93void G4SingleParticleSource::SetVerbosity(int vL)
94{
95  verbosityLevel = vL;
96  posGenerator->SetVerbosity(vL);
97  angGenerator->SetVerbosity(vL);
98  eneGenerator->SetVerbosity(vL);
99  G4cout << "Verbosity Set to: " << verbosityLevel << G4endl;
100}
101
102void G4SingleParticleSource::SetParticleDefinition
103  (G4ParticleDefinition* aParticleDefinition)
104{
105  particle_definition = aParticleDefinition;
106  particle_charge = particle_definition->GetPDGCharge();
107}
108
109void G4SingleParticleSource::GeneratePrimaryVertex(G4Event *evt)
110{
111  if(particle_definition==NULL) return;
112
113  if(verbosityLevel > 1)
114    G4cout << " NumberOfParticlesToBeGenerated: "<<NumberOfParticlesToBeGenerated << G4endl;
115
116  // Position stuff
117  particle_position = posGenerator->GenerateOne();
118
119  // create a new vertex
120  G4PrimaryVertex* vertex = new G4PrimaryVertex(particle_position,particle_time);
121
122  for( G4int i=0; i<NumberOfParticlesToBeGenerated; i++ ) {
123    // Angular stuff
124    particle_momentum_direction = angGenerator->GenerateOne();
125    // Energy stuff
126    particle_energy = eneGenerator->GenerateOne(particle_definition);
127   
128    if(verbosityLevel >= 2)
129      G4cout << "Creating primaries and assigning to vertex" << G4endl;
130    // create new primaries and set them to the vertex
131    G4double mass =  particle_definition->GetPDGMass();
132    G4double energy = particle_energy + mass;
133    G4double pmom = std::sqrt(energy*energy-mass*mass);
134    G4double px = pmom*particle_momentum_direction.x();
135    G4double py = pmom*particle_momentum_direction.y();
136    G4double pz = pmom*particle_momentum_direction.z();
137
138    if(verbosityLevel > 1){
139      G4cout << "Particle name: "<<particle_definition->GetParticleName() << G4endl; 
140      G4cout << "       Energy: "<<particle_energy << G4endl;
141      G4cout << "     Position: "<<particle_position<< G4endl; 
142      G4cout << "    Direction: "<<particle_momentum_direction << G4endl;
143    }
144    G4PrimaryParticle* particle =
145      new G4PrimaryParticle(particle_definition,px,py,pz);
146    particle->SetMass( mass );
147    particle->SetCharge( particle_charge );
148    particle->SetPolarization(particle_polarization.x(),
149                              particle_polarization.y(),
150                              particle_polarization.z());
151    // Set bweight equal to the multiple of all non-zero weights
152    particle_weight = biasRndm->GetBiasWeight();
153    // pass it to primary particle
154    particle->SetWeight(particle_weight);
155
156    vertex->SetPrimary( particle );
157     
158  }
159  // now pass the weight to the primary vertex. CANNOT be used here!
160  //  vertex->SetWeight(particle_weight);
161  evt->AddPrimaryVertex( vertex );
162  if(verbosityLevel > 1)
163    G4cout << " Primary Vetex generated !"<< G4endl;   
164}
165
166
167
168
169
170
171
172
173
Note: See TracBrowser for help on using the repository browser.