Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointhMultipleScattering.cc
Go to the documentation of this file.
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$
27//
28
29//
30// GEANT4 Class file
31//
32// File name: G4AdjointhMultipleScattering
33//
34// Author: Desorgher Laurent
35//
36// Creation date: 03.06.2009 cloned from G4hMultipleScattering by U.Laszlo with slight modification for adjoint_ion.
37//
38// -----------------------------------------------------------------------------
39//
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
44#include "G4SystemOfUnits.hh"
45#include "G4UrbanMscModel90.hh"
46#include "G4MscStepLimitType.hh"
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49
50using namespace std;
51
53 : G4VMultipleScattering(processName)
54{
55 isInitialized = false;
56 isIon = false;
58
59 dtrl=0.;
60 lambdalimit=0.;
61 samplez=0.;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
67{}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
72{
73 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77
79{
80 // Modification of parameters between runs
81 if(isInitialized) {
82 if (p->GetParticleType() != "adjoint_nucleus" && p->GetPDGMass() < GeV) {
85 mscUrban->SetSkin(Skin());
86 mscUrban->SetRangeFactor(RangeFactor());
87 mscUrban->SetGeomFactor(GeomFactor());
88 }
89 return;
90 }
91
92 // defaults for ions, which cannot be overwritten
93 if (p->GetParticleType() == "adjoint_nucleus" || p->GetPDGMass() > GeV) {
96 //SetBuildLambdaTable(false);
97 if(p->GetParticleType() == "adjoint_nucleus") isIon = true;
98 }
99
100 // initialisation of parameters
101 G4String part_name = p->GetParticleName();
102 mscUrban = new G4UrbanMscModel90();
103
104 mscUrban->SetStepLimitType(StepLimitType());
106 mscUrban->SetSkin(Skin());
107 mscUrban->SetRangeFactor(RangeFactor());
108 mscUrban->SetGeomFactor(GeomFactor());
109
110 AddEmModel(1,mscUrban);
111 isInitialized = true;
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
117{
118 G4cout << " RangeFactor= " << RangeFactor()
119 << ", step limit type: " << StepLimitType()
120 << ", lateralDisplacement: " << LateralDisplasmentFlag()
121 << ", skin= " << Skin()
122 // << ", geomFactor= " << GeomFactor()
123 << G4endl;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127
128/*G4double G4AdjointhMultipleScattering::AlongStepGetPhysicalInteractionLength(
129 const G4Track& track,
130 double,
131 G4double currentMinimalStep,
132 G4double& currentSafety,
133 G4GPILSelection* selection)
134{
135 // get Step limit proposed by the process
136 valueGPILSelectionMSC = NotCandidateForSelection;
137
138 G4double escaled = track.GetKineticEnergy();
139 if(isIon) escaled *= track.GetDynamicParticle()->GetMass()/proton_mass_c2;
140
141 G4double steplength = GetMscContinuousStepLimit(track,
142 escaled,
143 currentMinimalStep,
144 currentSafety);
145 // G4cout << "StepLimit= " << steplength << G4endl;
146 // set return value for G4GPILSelection
147 *selection = valueGPILSelectionMSC;
148 return steplength;
149}
150*/
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fMinimal
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void InitialiseProcess(const G4ParticleDefinition *)
G4bool IsApplicable(const G4ParticleDefinition &p)
G4AdjointhMultipleScattering(const G4String &processName="msc")
const G4String & GetParticleType() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void SetSkin(G4double)
Definition: G4VMscModel.hh:203
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:210
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:196
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:217
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:224
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetStepLimitType(G4MscStepLimitType val)
void SetLateralDisplasmentFlag(G4bool val)
G4bool LateralDisplasmentFlag() const
G4MscStepLimitType StepLimitType() const