source: trunk/source/processes/electromagnetic/adjoint/src/G4AdjointhMultipleScattering.cc@ 1259

Last change on this file since 1259 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 5.6 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: G4AdjointhMultipleScattering.cc,v 1.2 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29
30//
31// GEANT4 Class file
32//
33// File name: G4AdjointhMultipleScattering
34//
35// Author: Desorgher Laurent
36//
37// Creation date: 03.06.2009 cloned from G4hMultipleScattering by U.Laszlo with slight modification for adjoint_ion.
38//
39// -----------------------------------------------------------------------------
40//
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44#include "G4AdjointhMultipleScattering.hh"
45#include "G4UrbanMscModel.hh"
46#include "G4UrbanMscModel90.hh"
47#include "G4MscStepLimitType.hh"
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
51using namespace std;
52
53G4AdjointhMultipleScattering::G4AdjointhMultipleScattering(const G4String& processName)
54 : G4VMultipleScattering(processName)
55{
56 isInitialized = false;
57 isIon = false;
58 SetStepLimitType(fMinimal);
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
63G4AdjointhMultipleScattering::~G4AdjointhMultipleScattering()
64{}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68G4bool G4AdjointhMultipleScattering::IsApplicable (const G4ParticleDefinition& p)
69{
70 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
75void G4AdjointhMultipleScattering::InitialiseProcess(const G4ParticleDefinition* p)
76{
77 // Modification of parameters between runs
78 if(isInitialized) {
79 if (p->GetParticleType() != "adjoint_nucleus" && p->GetPDGMass() < GeV) {
80 mscUrban->SetStepLimitType(StepLimitType());
81 mscUrban->SetLateralDisplasmentFlag(LateralDisplasmentFlag());
82 mscUrban->SetSkin(Skin());
83 mscUrban->SetRangeFactor(RangeFactor());
84 mscUrban->SetGeomFactor(GeomFactor());
85 }
86 return;
87 }
88
89 // defaults for ions, which cannot be overwritten
90 if (p->GetParticleType() == "adjoint_nucleus" || p->GetPDGMass() > GeV) {
91 SetStepLimitType(fMinimal);
92 SetLateralDisplasmentFlag(false);
93 SetBuildLambdaTable(false);
94 if(p->GetParticleType() == "adjoint_nucleus") isIon = true;
95 }
96
97 // initialisation of parameters
98 G4String part_name = p->GetParticleName();
99 mscUrban = new G4UrbanMscModel90();
100
101 mscUrban->SetStepLimitType(StepLimitType());
102 mscUrban->SetLateralDisplasmentFlag(LateralDisplasmentFlag());
103 mscUrban->SetSkin(Skin());
104 mscUrban->SetRangeFactor(RangeFactor());
105 mscUrban->SetGeomFactor(GeomFactor());
106
107 AddEmModel(1,mscUrban);
108 isInitialized = true;
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
113void G4AdjointhMultipleScattering::PrintInfo()
114{
115 G4cout << " RangeFactor= " << RangeFactor()
116 << ", step limit type: " << StepLimitType()
117 << ", lateralDisplacement: " << LateralDisplasmentFlag()
118 << ", skin= " << Skin()
119 // << ", geomFactor= " << GeomFactor()
120 << G4endl;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
125/*G4double G4AdjointhMultipleScattering::AlongStepGetPhysicalInteractionLength(
126 const G4Track& track,
127 double,
128 G4double currentMinimalStep,
129 G4double& currentSafety,
130 G4GPILSelection* selection)
131{
132 // get Step limit proposed by the process
133 valueGPILSelectionMSC = NotCandidateForSelection;
134
135 G4double escaled = track.GetKineticEnergy();
136 if(isIon) escaled *= track.GetDynamicParticle()->GetMass()/proton_mass_c2;
137
138 G4double steplength = GetMscContinuousStepLimit(track,
139 escaled,
140 currentMinimalStep,
141 currentSafety);
142 // G4cout << "StepLimit= " << steplength << G4endl;
143 // set return value for G4GPILSelection
144 *selection = valueGPILSelectionMSC;
145 return steplength;
146}
147*/
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
Note: See TracBrowser for help on using the repository browser.