source: trunk/source/processes/electromagnetic/standard/src/G4hMultipleScattering.cc@ 899

Last change on this file since 899 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 4.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: G4hMultipleScattering.cc,v 1.7 2007/12/07 17:35:52 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
28//
29// -----------------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33// File name: G4hMultipleScattering
34//
35// Author: Laszlo Urban
36//
37// Creation date: 24.10.2006 cloned from G4MultipleScattering
38//
39// Modified:
40// 12-02-07 skin can be changed via UI command (VI)
41// 20.03.07 Remove local parameter skin, set facgeom=0.1(V.Ivanchenko)
42//
43// -----------------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
48#include "G4hMultipleScattering.hh"
49#include "G4UrbanMscModel.hh"
50#include "G4UrbanMscModel90.hh"
51#include "G4MscStepLimitType.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
55using namespace std;
56
57G4hMultipleScattering::G4hMultipleScattering(const G4String& processName)
58 : G4VMultipleScattering(processName)
59{
60 dtrl = 0.05;
61 lambdalimit = 1.*mm;
62
63 samplez = false ;
64 isInitialized = false;
65
66 SetLateralDisplasmentFlag(true);
67 SetSkin(0.0);
68 SetRangeFactor(0.2);
69 SetGeomFactor(0.1);
70 SetStepLimitType(fMinimal);
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
75G4hMultipleScattering::~G4hMultipleScattering()
76{}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
80G4bool G4hMultipleScattering::IsApplicable (const G4ParticleDefinition& p)
81{
82 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
87void G4hMultipleScattering::InitialiseProcess(const G4ParticleDefinition* p)
88{
89 // Modification of parameters between runs
90 if(isInitialized) {
91 if (p->GetParticleType() != "nucleus") {
92 mscUrban->SetStepLimitType(StepLimitType());
93 mscUrban->SetLateralDisplasmentFlag(LateralDisplasmentFlag());
94 mscUrban->SetSkin(Skin());
95 mscUrban->SetRangeFactor(RangeFactor());
96 mscUrban->SetGeomFactor(GeomFactor());
97 }
98 return;
99 }
100
101 // initialisation of parameters
102 G4String part_name = p->GetParticleName();
103 mscUrban = new G4UrbanMscModel90(RangeFactor(),dtrl,lambdalimit,
104 GeomFactor(),Skin(),
105 samplez,StepLimitType());
106 mscUrban->SetLateralDisplasmentFlag(LateralDisplasmentFlag());
107
108 if (p->GetParticleType() == "nucleus") {
109 mscUrban->SetStepLimitType(fMinimal);
110 SetLateralDisplasmentFlag(false);
111 SetBuildLambdaTable(false);
112 SetSkin(0.0);
113 SetRangeFactor(0.2);
114 }
115 AddEmModel(1,mscUrban);
116 isInitialized = true;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
121void G4hMultipleScattering::PrintInfo()
122{
123 G4cout << " Boundary/stepping algorithm is active with RangeFactor= "
124 << RangeFactor()
125 << " Step limit type " << StepLimitType()
126 << G4endl;
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130
Note: See TracBrowser for help on using the repository browser.