source: trunk/source/processes/electromagnetic/standard/src/G4MultipleScattering71.cc @ 825

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

import all except CVS

File size: 7.7 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: G4MultipleScattering71.cc,v 1.4 2006/10/16 15:26:49 vnivanch Exp $
27// GEANT4 tag $Name:  $
28//
29// -----------------------------------------------------------------------------
30// 16/05/01 value of cparm changed , L.Urban
31// 18/05/01 V.Ivanchenko Clean up against Linux ANSI compilation
32// 07/08/01 new methods Store/Retrieve PhysicsTable (mma)
33// 23-08-01 new angle and z distribution,energy dependence reduced,
34//          Store,Retrieve methods commented out temporarily, L.Urban
35// 27-08-01 in BuildPhysicsTable:aParticleType.GetParticleName()=="mu+" (mma)
36// 28-08-01 GetContinuousStepLimit and AlongStepDoIt moved from .icc file (mma)
37// 03-09-01 value of data member factlim changed, L.Urban
38// 10-09-01 small change in GetContinuousStepLimit, L.Urban
39// 11-09-01 G4MultipleScatteringx put as default G4MultipleScattering
40//          store/retrieve physics table reactivated (mma)
41// 13-09-01 corr. in ComputeTransportCrossSection, L.Urban
42// 14-09-01 protection in GetContinuousStepLimit, L.Urban
43// 17-09-01 migration of Materials to pure STL (mma)
44// 27-09-01 value of data member factlim changed, L.Urban
45// 31-10-01 big fixed in PostStepDoIt,L.Urban
46// 17-04-02 NEW angle distribution + boundary algorithm modified, L.Urban
47// 22-04-02 boundary algorithm modified -> important improvement in timing (L.Urban)
48// 24-04-02 some minor changes in boundary algorithm, L.Urban
49// 06-05-02 bug fixed in GetContinuousStepLimit, L.Urban
50// 24-05-02 changes in angle distribution and boundary algorithm, L.Urban
51// 11-06-02 bug fixed in ComputeTransportCrossSection, L.Urban
52// 12-08-02 bug fixed in PostStepDoIt (lateral displacement), L.Urban
53// 15-08-02 new angle distribution, L.Urban
54// 26-09-02 angle distribution + boundary algorithm modified, L.Urban
55// 15-10-02 temporary fix for proton scattering
56// 30-10-02 modified angle distribution,mods in boundary algorithm,
57//          changes in data members, L.Urban
58// 11-12-02 precision problem in ComputeTransportCrossSection
59//          for small Tkin/for heavy particles cured from L.Urban
60// 20-01-03 Migrade to cut per region (V.Ivanchenko)
61// 05-02-03 changes in data members, new sampling for geom.
62//          path length, step dependence reduced with new
63//          method
64// 28-03-03 Move to model design (V.Ivanchenko)
65// 08-08-03 STD substitute standard  (V.Ivanchenko)
66// 23-04-04 value of data member dtrl changed from 0.15 to 0.05 (L.Urban)
67// 17-08-04 name of facxsi changed to factail (L.Urban)
68// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
69// 07-02-05 correction in order to have a working Setsamplez function (L.Urban)
70// 15-04-05 optimize internal interface (V.Ivanchenko)
71// 03-10-05 Process is freezed with the name 71 (V.Ivanchenko)
72// -----------------------------------------------------------------------------
73//
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
77#include "G4MultipleScattering71.hh"
78#include "G4MscModel71.hh"
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
82using namespace std;
83
84G4MultipleScattering71::G4MultipleScattering71(const G4String& processName)
85  : G4VMultipleScattering(processName),
86    totBins(120),
87    facrange(0.199),
88    dtrl(0.05),
89    NuclCorrPar (0.0615),
90    FactPar(0.40),
91    factail(1.0),
92    cf(1.001),
93    stepnolastmsc(-1000000),
94    nsmallstep(5),
95    samplez(true),
96    boundary(true),
97    isInitialized(false)
98{
99  lowKineticEnergy = 0.1*keV;
100  highKineticEnergy= 100.*TeV;
101
102  tlimit           = 1.e10*mm;
103  tlimitmin        = 1.e-7*mm;
104
105  SetBinning(totBins);
106  SetMinKinEnergy(lowKineticEnergy);
107  SetMaxKinEnergy(highKineticEnergy);
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111
112G4MultipleScattering71::~G4MultipleScattering71()
113{}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
117void G4MultipleScattering71::InitialiseProcess(const G4ParticleDefinition* particle)
118{
119  if(isInitialized) return;
120
121  if (particle->GetParticleType() == "nucleus") {
122    boundary = false;
123    SetLateralDisplasmentFlag(false);
124    SetBuildLambdaTable(false);
125    Setsamplez(false) ;
126  } else {
127    SetLateralDisplasmentFlag(true);
128    SetBuildLambdaTable(true);
129  }
130  G4MscModel71* em = new G4MscModel71(dtrl,NuclCorrPar,FactPar,factail,samplez);
131  em->SetLateralDisplasmentFlag(LateralDisplasmentFlag());
132  em->SetLowEnergyLimit(lowKineticEnergy);
133  em->SetHighEnergyLimit(highKineticEnergy);
134  AddEmModel(1, em);
135  isInitialized = true;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140G4double G4MultipleScattering71::TruePathLengthLimit(const G4Track&  track,
141                                                            G4double& lambda,
142                                                            G4double  currentMinimalStep)
143{
144
145  G4double tPathLength = currentMinimalStep;
146
147  // special treatment near boundaries ?
148  if (boundary) {
149
150    G4int stepno = track.GetCurrentStepNumber() ;
151    // first step
152    if (stepno == 1) {
153      stepnolastmsc = -1000000 ;
154      tlimit = 1.e10;
155    } else if (stepno > 1) {
156
157      if (track.GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) {
158
159        stepnolastmsc = stepno;
160        //  if : diff.treatment for small/not small Z
161        if (range > lambda) tlimit = facrange*range;
162        else                tlimit = facrange*lambda;
163
164        if(tlimit < tlimitmin) tlimit = tlimitmin;
165        if(tPathLength > tlimit) tPathLength = tlimit;
166
167      } else if (stepno > stepnolastmsc && stepno - stepnolastmsc < nsmallstep
168              && tPathLength > tlimit) {
169        tlimit *= cf;
170        tPathLength = tlimit;
171      }
172    }
173  }
174
175  return tPathLength;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
180void G4MultipleScattering71::PrintInfo()
181{
182  if(boundary) {
183    G4cout << "      Boundary algorithm is active with facrange= "
184           << facrange
185           << G4endl;
186  }
187  G4cout << "        WARNING: This process is obsolete and will be soon removed" 
188         << G4endl;
189 
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193
Note: See TracBrowser for help on using the repository browser.