source: trunk/source/processes/electromagnetic/utils/src/G4VMultipleScattering.cc@ 968

Last change on this file since 968 was 961, checked in by garnier, 17 years ago

update processes

File size: 13.8 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// $Id: G4VMultipleScattering.cc,v 1.60 2008/11/20 20:32:40 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4VMultipleScattering
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 25.03.2003
39//
40// Modifications:
41//
42// 13.04.03 Change printout (V.Ivanchenko)
43// 04-06-03 Fix compilation warnings (V.Ivanchenko)
44// 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
45// 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
46// 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
47// 01-03-04 SampleCosineTheta signature changed
48// 22-04-04 SampleCosineTheta signature changed back to original
49// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
50// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
51// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
52// 15-04-05 optimize internal interface (V.Ivanchenko)
53// 15-04-05 remove boundary flag (V.Ivanchenko)
54// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
55// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
56// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
57// 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
58//
59// Class Description:
60//
61// It is the generic process of multiple scattering it includes common
62// part of calculations for all charged particles
63
64// -------------------------------------------------------------------
65//
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
69#include "G4VMultipleScattering.hh"
70#include "G4LossTableManager.hh"
71#include "G4Step.hh"
72#include "G4ParticleDefinition.hh"
73#include "G4VEmFluctuationModel.hh"
74#include "G4DataVector.hh"
75#include "G4PhysicsTable.hh"
76#include "G4PhysicsVector.hh"
77#include "G4PhysicsLogVector.hh"
78#include "G4UnitsTable.hh"
79#include "G4ProductionCutsTable.hh"
80#include "G4Region.hh"
81#include "G4RegionStore.hh"
82#include "G4PhysicsTableHelper.hh"
83#include "G4GenericIon.hh"
84#include "G4Electron.hh"
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
88G4VMultipleScattering::G4VMultipleScattering(const G4String& name,
89 G4ProcessType type):
90 G4VContinuousDiscreteProcess(name, type),
91 buildLambdaTable(true),
92 theLambdaTable(0),
93 firstParticle(0),
94 stepLimit(fUseSafety),
95 skin(3.0),
96 facrange(0.02),
97 facgeom(2.5),
98 latDisplasment(true),
99 currentParticle(0),
100 currentCouple(0)
101{
102 SetVerboseLevel(1);
103 SetProcessSubType(fMultipleScattering);
104
105 // Size of tables assuming spline
106 minKinEnergy = 0.1*keV;
107 maxKinEnergy = 100.0*TeV;
108 nBins = 84;
109
110 // default limit on polar angle
111 polarAngleLimit = 0.0;
112
113 pParticleChange = &fParticleChange;
114
115 modelManager = new G4EmModelManager();
116 (G4LossTableManager::Instance())->Register(this);
117
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
122G4VMultipleScattering::~G4VMultipleScattering()
123{
124 if(1 < verboseLevel) {
125 G4cout << "G4VMultipleScattering destruct " << GetProcessName()
126 << G4endl;
127 }
128 delete modelManager;
129 if (theLambdaTable) {
130 theLambdaTable->clearAndDestroy();
131 delete theLambdaTable;
132 }
133 (G4LossTableManager::Instance())->DeRegister(this);
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138void G4VMultipleScattering::BuildPhysicsTable(const G4ParticleDefinition& part)
139{
140 G4String num = part.GetParticleName();
141 if(1 < verboseLevel) {
142 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
143 << GetProcessName()
144 << " and particle " << num
145 << G4endl;
146 }
147
148 if (buildLambdaTable && firstParticle == &part) {
149
150 const G4ProductionCutsTable* theCoupleTable=
151 G4ProductionCutsTable::GetProductionCutsTable();
152 size_t numOfCouples = theCoupleTable->GetTableSize();
153
154 for (size_t i=0; i<numOfCouples; i++) {
155
156 if (theLambdaTable->GetFlag(i)) {
157 // create physics vector and fill it
158 const G4MaterialCutsCouple* couple =
159 theCoupleTable->GetMaterialCutsCouple(i);
160 G4PhysicsVector* aVector = PhysicsVector(couple);
161 modelManager->FillLambdaVector(aVector, couple, false);
162 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector);
163 }
164 }
165
166 if(1 < verboseLevel) {
167 G4cout << "Lambda table is built for "
168 << num
169 << G4endl;
170 }
171 }
172 if(verboseLevel>0 && ( num == "e-" || num == "mu+" ||
173 num == "proton" || num == "pi-" ||
174 num == "GenericIon")) {
175 PrintInfoDefinition();
176 if(2 < verboseLevel && theLambdaTable) G4cout << *theLambdaTable << G4endl;
177 }
178
179 if(1 < verboseLevel) {
180 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
181 << GetProcessName()
182 << " and particle " << num
183 << G4endl;
184 }
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189void G4VMultipleScattering::PreparePhysicsTable(const G4ParticleDefinition& part)
190{
191 if (!firstParticle) {
192 currentCouple = 0;
193 if(part.GetParticleType() == "nucleus" &&
194 part.GetParticleSubType() == "generic") {
195 firstParticle = G4GenericIon::GenericIon();
196 } else {
197 firstParticle = &part;
198 }
199
200 currentParticle = &part;
201 }
202
203 if(1 < verboseLevel) {
204 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
205 << GetProcessName()
206 << " and particle " << part.GetParticleName()
207 << " local particle " << firstParticle->GetParticleName()
208 << G4endl;
209 }
210
211 if(firstParticle == &part) {
212
213 InitialiseProcess(firstParticle);
214 if(buildLambdaTable)
215 theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
216 const G4DataVector* theCuts =
217 modelManager->Initialise(firstParticle,
218 G4Electron::Electron(),
219 10.0, verboseLevel);
220
221 if(2 < verboseLevel) G4cout << theCuts << G4endl;
222
223 }
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
228void G4VMultipleScattering::PrintInfoDefinition()
229{
230 if (0 < verboseLevel) {
231 G4cout << G4endl << GetProcessName()
232 << ": for " << firstParticle->GetParticleName()
233 << " SubType= " << GetProcessSubType()
234 << G4endl;
235 if (theLambdaTable) {
236 G4cout << " Lambda tables from "
237 << G4BestUnit(MinKinEnergy(),"Energy")
238 << " to "
239 << G4BestUnit(MaxKinEnergy(),"Energy")
240 << " in " << nBins << " bins, spline: "
241 << (G4LossTableManager::Instance())->SplineFlag()
242 << G4endl;
243 }
244 PrintInfo();
245 modelManager->DumpModelList(verboseLevel);
246 if (2 < verboseLevel) {
247 G4cout << "LambdaTable address= " << theLambdaTable << G4endl;
248 if(theLambdaTable) G4cout << (*theLambdaTable) << G4endl;
249 }
250 }
251}
252
253//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
254
255G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(
256 const G4Track& track,
257 G4double,
258 G4double currentMinimalStep,
259 G4double& currentSafety,
260 G4GPILSelection* selection)
261{
262 // get Step limit proposed by the process
263 valueGPILSelectionMSC = NotCandidateForSelection;
264 G4double steplength = GetMscContinuousStepLimit(track,
265 track.GetKineticEnergy(),
266 currentMinimalStep,
267 currentSafety);
268 // G4cout << "StepLimit= " << steplength << G4endl;
269 // set return value for G4GPILSelection
270 *selection = valueGPILSelectionMSC;
271 return steplength;
272}
273
274//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
275
276G4double G4VMultipleScattering::PostStepGetPhysicalInteractionLength(
277 const G4Track&, G4double, G4ForceCondition* condition)
278{
279 *condition = Forced;
280 return DBL_MAX;
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
284
285G4double G4VMultipleScattering::GetContinuousStepLimit(
286 const G4Track& track,
287 G4double previousStepSize,
288 G4double currentMinimalStep,
289 G4double& currentSafety)
290{
291 return GetMscContinuousStepLimit(track,previousStepSize,currentMinimalStep,
292 currentSafety);
293}
294
295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
296
297G4double G4VMultipleScattering::GetMeanFreePath(
298 const G4Track&, G4double, G4ForceCondition* condition)
299{
300 *condition = Forced;
301 return DBL_MAX;
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
305
306G4VParticleChange* G4VMultipleScattering::AlongStepDoIt(const G4Track&,
307 const G4Step& step)
308{
309 fParticleChange.ProposeTrueStepLength(
310 currentModel->ComputeTrueStepLength(step.GetStepLength()));
311 return &fParticleChange;
312}
313
314//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
315
316G4VParticleChange* G4VMultipleScattering::PostStepDoIt(const G4Track& track,
317 const G4Step& step)
318{
319 fParticleChange.Initialize(track);
320 currentModel->SampleScattering(track.GetDynamicParticle(),
321 step.GetPostStepPoint()->GetSafety());
322 return &fParticleChange;
323}
324
325//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
326
327G4PhysicsVector* G4VMultipleScattering::PhysicsVector(const G4MaterialCutsCouple* couple)
328{
329 G4int nbins = 3;
330 if( couple->IsUsed() ) nbins = nBins;
331 G4PhysicsVector* v = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nbins);
332 v->SetSpline((G4LossTableManager::Instance())->SplineFlag());
333 return v;
334}
335
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
337
338G4bool G4VMultipleScattering::StorePhysicsTable(const G4ParticleDefinition* part,
339 const G4String& directory,
340 G4bool ascii)
341{
342 G4bool yes = true;
343 if ( theLambdaTable && part == firstParticle) {
344 const G4String name = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
345 G4bool yes = theLambdaTable->StorePhysicsTable(name,ascii);
346
347 if ( yes ) {
348 if ( verboseLevel>0 ) {
349 G4cout << "Physics table are stored for " << part->GetParticleName()
350 << " and process " << GetProcessName()
351 << " in the directory <" << directory
352 << "> " << G4endl;
353 }
354 } else {
355 G4cout << "Fail to store Physics Table for " << part->GetParticleName()
356 << " and process " << GetProcessName()
357 << " in the directory <" << directory
358 << "> " << G4endl;
359 }
360 }
361 return yes;
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
365
366G4bool G4VMultipleScattering::RetrievePhysicsTable(const G4ParticleDefinition* part,
367 const G4String& directory,
368 G4bool ascii)
369{
370 if(0 < verboseLevel) {
371// G4cout << "========================================================" << G4endl;
372 G4cout << "G4VMultipleScattering::RetrievePhysicsTable() for "
373 << part->GetParticleName() << " and process "
374 << GetProcessName() << G4endl;
375 }
376 G4bool yes = true;
377
378 if(!buildLambdaTable || firstParticle != part) return yes;
379
380 const G4String particleName = part->GetParticleName();
381
382 G4String filename = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
383 yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTable,filename,ascii);
384 if ( yes ) {
385 if (0 < verboseLevel) {
386 G4cout << "Lambda table for " << part->GetParticleName() << " is retrieved from <"
387 << filename << ">"
388 << G4endl;
389 }
390 if((G4LossTableManager::Instance())->SplineFlag()) {
391 size_t n = theLambdaTable->length();
392 for(size_t i=0; i<n; i++) {(* theLambdaTable)[i]->SetSpline(true);}
393 }
394 } else {
395 if (1 < verboseLevel) {
396 G4cout << "Lambda table for " << part->GetParticleName() << " in file <"
397 << filename << "> is not exist"
398 << G4endl;
399 }
400 }
401 return yes;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
405
Note: See TracBrowser for help on using the repository browser.