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

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

import all except CVS

File size: 13.5 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.47 2007/11/09 11:35:54 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-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//
58// Class Description:
59//
60// It is the generic process of multiple scattering it includes common
61// part of calculations for all charged particles
62
63// -------------------------------------------------------------------
64//
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
68#include "G4VMultipleScattering.hh"
69#include "G4LossTableManager.hh"
70#include "G4Step.hh"
71#include "G4ParticleDefinition.hh"
72#include "G4VEmFluctuationModel.hh"
73#include "G4DataVector.hh"
74#include "G4PhysicsTable.hh"
75#include "G4PhysicsVector.hh"
76#include "G4PhysicsLogVector.hh"
77#include "G4UnitsTable.hh"
78#include "G4ProductionCutsTable.hh"
79#include "G4Region.hh"
80#include "G4RegionStore.hh"
81#include "G4PhysicsTableHelper.hh"
82#include "G4GenericIon.hh"
83#include "G4Electron.hh"
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
87G4VMultipleScattering::G4VMultipleScattering(const G4String& name, G4ProcessType type):
88                 G4VContinuousDiscreteProcess(name, type),
89  theLambdaTable(0),
90  firstParticle(0),
91  currentParticle(0),
92  currentCouple(0),
93  nBins(120),
94  stepLimit(fUseSafety),
95  skin(0.0),
96  facrange(0.02),
97  facgeom(2.5),
98  latDisplasment(true),
99  buildLambdaTable(true)
100{
101  minKinEnergy = 0.1*keV;
102  maxKinEnergy = 100.0*TeV;
103  SetVerboseLevel(1);
104
105  pParticleChange = &fParticleChange;
106
107  modelManager = new G4EmModelManager();
108  (G4LossTableManager::Instance())->Register(this);
109
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
114G4VMultipleScattering::~G4VMultipleScattering()
115{
116  if(1 < verboseLevel) 
117    G4cout << "G4VMultipleScattering destruct " << GetProcessName() 
118           << G4endl;
119  delete modelManager;
120  if (theLambdaTable) {
121    theLambdaTable->clearAndDestroy();
122    delete theLambdaTable;
123  }
124  (G4LossTableManager::Instance())->DeRegister(this);
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
129void G4VMultipleScattering::BuildPhysicsTable(const G4ParticleDefinition& part)
130{
131  G4String num = part.GetParticleName();
132  if(1 < verboseLevel) {
133    //    G4cout << "========================================================" << G4endl;
134    G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
135           << GetProcessName()
136           << " and particle " << num
137           << G4endl;
138  }
139
140  if (buildLambdaTable && firstParticle == &part) {
141
142    const G4ProductionCutsTable* theCoupleTable=
143          G4ProductionCutsTable::GetProductionCutsTable();
144    size_t numOfCouples = theCoupleTable->GetTableSize();
145
146    for (size_t i=0; i<numOfCouples; i++) {
147
148      if (theLambdaTable->GetFlag(i)) {
149        // create physics vector and fill it
150        const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
151        G4PhysicsVector* aVector = PhysicsVector(couple);
152        modelManager->FillLambdaVector(aVector, couple, false);
153        G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector);
154      }
155    }
156
157    if(1 < verboseLevel) {
158      G4cout << "Lambda table is built for "
159             << num
160             << G4endl;
161    }
162  }
163  if(verboseLevel>0 && ( num == "e-" || num == "mu+" || 
164                         num == "proton" || num == "pi-" || num == "GenericIon")) {
165    PrintInfoDefinition();
166    if(2 < verboseLevel && theLambdaTable) G4cout << *theLambdaTable << G4endl;
167  }
168
169  if(1 < verboseLevel) {
170    G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
171           << GetProcessName()
172           << " and particle " << num
173           << G4endl;
174  }
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
179void G4VMultipleScattering::PreparePhysicsTable(const G4ParticleDefinition& part)
180{
181  if (!firstParticle) {
182    currentCouple = 0;
183    if(part.GetParticleType() == "nucleus" && 
184       part.GetParticleSubType() == "generic") 
185         firstParticle = G4GenericIon::GenericIon();
186    else firstParticle = &part; 
187    currentParticle = &part;
188  }
189
190  if(1 < verboseLevel) {
191    //    G4cout << "========================================================" << G4endl;
192    G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
193           << GetProcessName()
194           << " and particle " << part.GetParticleName()
195           << " local particle " << firstParticle->GetParticleName()
196           << G4endl;
197  }
198
199  if(firstParticle == &part) {
200
201    InitialiseProcess(firstParticle);
202    if(buildLambdaTable)
203      theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
204    const G4DataVector* theCuts = 
205      modelManager->Initialise(firstParticle, 
206                               G4Electron::Electron(), 
207                               10.0, verboseLevel);
208
209    if(2 < verboseLevel) G4cout << theCuts << G4endl;
210
211  }
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215
216void G4VMultipleScattering::PrintInfoDefinition()
217{
218  if (0 < verboseLevel) {
219    G4cout << G4endl << GetProcessName() << ":  Model variant of multiple scattering "
220           << "for " << firstParticle->GetParticleName()
221           << G4endl;
222    if (theLambdaTable) {
223      G4cout << "      Lambda tables from "
224             << G4BestUnit(MinKinEnergy(),"Energy")
225             << " to "
226             << G4BestUnit(MaxKinEnergy(),"Energy")
227             << " in " << nBins << " bins."
228             << G4endl;
229    }
230    G4cout << "      LateralDisplacementFlag=  " << latDisplasment
231           << "   Skin= " << skin << G4endl; 
232    PrintInfo();
233    if (2 < verboseLevel) {
234      G4cout << "LambdaTable address= " << theLambdaTable << G4endl;
235      if(theLambdaTable) G4cout << (*theLambdaTable) << G4endl;
236    }
237  }
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
241
242G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(
243                             const G4Track& track,
244                             G4double previousStepSize,
245                             G4double currentMinimalStep,
246                             G4double& currentSafety,
247                             G4GPILSelection* selection)
248{
249  // get Step limit proposed by the process
250  valueGPILSelectionMSC = NotCandidateForSelection;
251  G4double steplength = GetMscContinuousStepLimit(track,previousStepSize,
252                                              currentMinimalStep,currentSafety);
253  // G4cout << "StepLimit= " << steplength << G4endl;
254  // set return value for G4GPILSelection
255  *selection = valueGPILSelectionMSC;
256  return  steplength;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260
261G4double G4VMultipleScattering::PostStepGetPhysicalInteractionLength(
262              const G4Track&, G4double, G4ForceCondition* condition)
263{
264  *condition = Forced;
265  return DBL_MAX;
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
269
270G4double G4VMultipleScattering::GetContinuousStepLimit(
271                                       const G4Track& track,
272                                       G4double previousStepSize,
273                                       G4double currentMinimalStep,
274                                       G4double& currentSafety)
275{
276  return GetMscContinuousStepLimit(track,previousStepSize,currentMinimalStep,
277                                   currentSafety);
278}
279
280//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
281
282G4double G4VMultipleScattering::GetMeanFreePath(
283              const G4Track&, G4double, G4ForceCondition* condition)
284{
285  *condition = Forced;
286  return DBL_MAX;
287}
288
289//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
290
291G4VParticleChange* G4VMultipleScattering::AlongStepDoIt(const G4Track&,
292                                                        const G4Step& step)
293{
294  fParticleChange.ProposeTrueStepLength(
295    currentModel->ComputeTrueStepLength(step.GetStepLength()));
296  return &fParticleChange;
297}
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300
301G4VParticleChange* G4VMultipleScattering::PostStepDoIt(const G4Track& track,
302                                                       const G4Step& step)
303{
304  fParticleChange.Initialize(track);
305  currentModel->SampleScattering(track.GetDynamicParticle(),
306                                 step.GetPostStepPoint()->GetSafety());
307  return &fParticleChange;
308}
309
310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311
312G4PhysicsVector* G4VMultipleScattering::PhysicsVector(const G4MaterialCutsCouple* couple)
313{
314  G4int nbins = 3;
315  if( couple->IsUsed() ) nbins = nBins;
316  G4PhysicsVector* v = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nbins);
317  return v;
318}
319
320//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
321
322G4bool G4VMultipleScattering::StorePhysicsTable(const G4ParticleDefinition* part,
323                                                const G4String& directory,
324                                                      G4bool ascii)
325{
326  G4bool yes = true;
327  if ( theLambdaTable && part == firstParticle) {
328    const G4String name = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
329    G4bool yes = theLambdaTable->StorePhysicsTable(name,ascii);
330
331    if ( yes ) {
332      if ( verboseLevel>0 ) {
333        G4cout << "Physics table are stored for " << part->GetParticleName()
334               << " and process " << GetProcessName()
335               << " in the directory <" << directory
336               << "> " << G4endl;
337      }
338    } else {
339      G4cout << "Fail to store Physics Table for " << part->GetParticleName()
340             << " and process " << GetProcessName()
341             << " in the directory <" << directory
342             << "> " << G4endl;
343    }
344  }
345  return yes;
346}
347
348//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
349
350G4bool G4VMultipleScattering::RetrievePhysicsTable(const G4ParticleDefinition* part,
351                                                   const G4String& directory,
352                                                         G4bool ascii)
353{
354  if(0 < verboseLevel) {
355//    G4cout << "========================================================" << G4endl;
356    G4cout << "G4VMultipleScattering::RetrievePhysicsTable() for "
357           << part->GetParticleName() << " and process "
358           << GetProcessName() << G4endl;
359  }
360  G4bool yes = true;
361
362  if(!buildLambdaTable || firstParticle != part) return yes;
363
364  const G4String particleName = part->GetParticleName();
365
366  G4String filename = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
367  yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTable,filename,ascii);
368  if ( yes ) {
369    if (0 < verboseLevel) {
370        G4cout << "Lambda table for " << part->GetParticleName() << " is retrieved from <"
371               << filename << ">"
372               << G4endl;
373    }
374  } else {
375    if (1 < verboseLevel) {
376        G4cout << "Lambda table for " << part->GetParticleName() << " in file <"
377               << filename << "> is not exist"
378               << G4endl;
379    }
380  }
381  return yes;
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
385
Note: See TracBrowser for help on using the repository browser.