source: trunk/source/processes/electromagnetic/lowenergy/test/processTest/src/G4ProcessTest.cc@ 1201

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

nvx fichiers dans CVS

File size: 6.3 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: G4ProcessTest.cc,v 1.12 2006/06/29 19:48:50 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30// Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// History:
33// -----------
34// 07 Oct 2001 MGP Created
35//
36// -------------------------------------------------------------------
37// Class description:
38// Test DoIt method of electromagnetic physics processes
39// Further documentation available from http://www.ge.infn.it/geant4/lowE/index.html
40
41// -------------------------------------------------------------------
42
43#include "globals.hh"
44#include "G4ProcessTest.hh"
45#include "G4VProcess.hh"
46#include "G4ProcessManager.hh"
47#include "G4Track.hh"
48#include "G4Step.hh"
49#include "G4VParticleChange.hh"
50#include "G4ParticleChange.hh"
51#include "G4ParticleDefinition.hh"
52#include "G4Gamma.hh"
53#include "G4Electron.hh"
54#include "G4Positron.hh"
55#include "G4eBremsstrahlung.hh"
56#include "G4eIonisation.hh"
57#include "G4eplusAnnihilation.hh"
58
59#include "G4ProcessTestAnalysis.hh"
60
61G4ProcessTest::G4ProcessTest() :
62 process(0), ioni(0), brem(0), eProcessManager(0), gProcessManager(0), def(0)
63{ }
64
65G4ProcessTest:: ~G4ProcessTest()
66{
67 delete process;
68 process = 0;
69 delete ioni;
70 ioni = 0;
71 delete brem;
72 brem = 0;
73 delete eProcessManager;
74 eProcessManager = 0;
75 delete gProcessManager;
76 gProcessManager = 0;
77 }
78
79void G4ProcessTest::buildTables(const G4String& type, G4bool isPolarised)
80{
81 G4ParticleDefinition* gamma = G4Gamma::GammaDefinition();
82 gamma->SetCuts(1e-3*mm);
83 G4ParticleDefinition* electron = G4Electron::ElectronDefinition();
84 electron->SetCuts(1e-3*mm);
85
86 def = createIncidentParticle();
87
88 process = createProcess();
89 G4cout << process->GetProcessName() << " created" << G4endl;
90
91 brem = createBremsstrahlung();
92 if (brem != 0) G4cout << brem->GetProcessName() << " created" << G4endl;
93
94 ioni = createElectronIonisation();
95 if (ioni != 0) G4cout << ioni->GetProcessName() << " created" << G4endl;
96
97 eProcessManager = new G4ProcessManager(electron);
98 electron->SetProcessManager(eProcessManager);
99
100 G4cout << "Now building physics tables, it will take a while..." << G4endl;
101
102 if (def == G4Gamma::GammaDefinition())
103 {
104 gProcessManager = new G4ProcessManager(def);
105 def->SetProcessManager(gProcessManager);
106 gProcessManager->AddProcess(process);
107 process->BuildPhysicsTable(*def);
108 }
109
110 if (def == G4Electron::ElectronDefinition())
111 {
112 eProcessManager->AddProcess(process);
113 process->BuildPhysicsTable(*def);
114 }
115
116 // Electron processes are always created; they are needed in photon
117 // processes for range tests
118 if (ioni != 0) eProcessManager->AddProcess(ioni);
119 if (brem != 0) eProcessManager->AddProcess(brem);
120 if (ioni != 0) ioni->BuildPhysicsTable(*electron);
121 if (brem != 0) brem->BuildPhysicsTable(*electron);
122
123 // Build (standard) physics tables for positrons, needed by GammaConversion
124
125 G4ParticleDefinition* positron = G4Positron::PositronDefinition();
126 positron->SetCuts(1e-3*mm);
127 G4VProcess* theeplusIonisation = new G4eIonisation();
128 G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
129 G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
130 G4ProcessManager* posProcessManager = new G4ProcessManager(positron);
131 positron->SetProcessManager(posProcessManager);
132 posProcessManager->AddProcess(theeplusIonisation);
133 posProcessManager->AddProcess(theeplusBremsstrahlung);
134 posProcessManager->AddProcess(theeplusAnnihilation);
135 theeplusIonisation->BuildPhysicsTable(*positron);
136 theeplusBremsstrahlung->BuildPhysicsTable(*positron);
137 theeplusAnnihilation->BuildPhysicsTable(*positron) ;
138}
139
140
141void G4ProcessTest::postStepTest(const G4Track& track,
142 const G4Step& step) const
143{
144 G4ParticleChange* particleChange = (G4ParticleChange*) process->PostStepDoIt(track, step);
145
146 G4ProcessTestAnalysis* analysis = G4ProcessTestAnalysis::getInstance();
147 analysis->analyseGeneral(track,particleChange);
148 analysis->analyseSecondaries(particleChange);
149
150 for (G4int i = 0; i < (particleChange->GetNumberOfSecondaries()); i++)
151 {
152 delete particleChange->GetSecondary(i);
153 }
154
155 particleChange->Clear();
156}
157
158void G4ProcessTest::alongStepTest(const G4Track& track,
159 const G4Step& step) const
160{
161 G4ParticleChange* particleChange = (G4ParticleChange*) process->PostStepDoIt(track, step);
162
163 G4ProcessTestAnalysis* analysis = G4ProcessTestAnalysis::getInstance();
164 analysis->analyseGeneral(track,particleChange);
165 analysis->analyseSecondaries(particleChange);
166
167 for (G4int i = 0; i < (particleChange->GetNumberOfSecondaries()); i++)
168 {
169 delete particleChange->GetSecondary(i);
170 }
171
172 particleChange->Clear();
173}
174
175G4ParticleDefinition* G4ProcessTest::createIncidentParticle()
176{
177 return G4Gamma::GammaDefinition();
178}
Note: See TracBrowser for help on using the repository browser.