OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
ScatteringPhysics.h
Go to the documentation of this file.
1//
2// Class ScatteringPhysics
3// Defines the physical processes of beam scattering
4// and energy loss by heavy charged particles
5//
6// Copyright (c) 2009 - 2022, Bi, Yang, Stachel, Adelmann
7// Paul Scherrer Institut, Villigen PSI, Switzerland
8// All rights reserved.
9//
10// This file is part of OPAL.
11//
12// OPAL is free software: you can redistribute it and/or modify
13// it under the terms of the GNU General Public License as published by
14// the Free Software Foundation, either version 3 of the License, or
15// (at your option) any later version.
16//
17// You should have received a copy of the GNU General Public License
18// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
19//
20#ifndef SCATTERINGPHYSICS_H
21#define SCATTERINGPHYSICS_H
22
24
26#include "Algorithms/Vektor.h"
27
28#include <gsl/gsl_rng.h>
29
30#include "Utility/IpplTimings.h"
31
32#include <memory>
33#include <string>
34#include <utility>
35#include <vector>
36
37template <class T, unsigned Dim>
38class PartBunchBase;
39
40class LossDataSink;
41class Inform;
42
43typedef struct { // struct for description of particle in material
44 int label; // the status of the particle (0 = in material / -1 = move back to bunch
45 unsigned localID; // not so unique identifier of the particle
46 Vector_t Rincol; // position
47 Vector_t Pincol; // momentum
48 long IDincol; // unique identifier of the particle inherited from the bunch
49 int Binincol; // bin number
50 double DTincol; // time step size
51 double Qincol; // charge
52 double Mincol; // mass
53 Vector_t Bfincol; // magnetic field
54 Vector_t Efincol; // electric field
55} PART;
56
57
59public:
60
61 ScatteringPhysics(const std::string& name,
62 ElementBase* element,
63 std::string& mat,
64 bool enableRutherford,
65 double lowEnergyThr);
67
68 virtual void apply(PartBunchBase<double, 3>* bunch,
69 const std::pair<Vector_t, double>& boundingSphere) override;
70
71 virtual const std::string getType() const override;
72 virtual void print(Inform& os) override;
73 virtual bool stillActive() override;
74
75 virtual double getTime() override;
76 virtual std::string getName() override;
77 virtual size_t getParticlesInMat() override;
78 virtual unsigned getRediffused() override;
79 virtual unsigned int getNumEntered() override;
80
82
84 Vector_t& P,
85 const double deltat,
86 bool includeFluctuations = true) const override;
87
88private:
89
92 Vector_t& P,
93 double dt);
94
95 void applyRotation(Vector_t& P,
96 Vector_t& R,
97 double xplane,
98 double thetacou);
99 void applyRandomRotation(Vector_t& P, double theta0);
100
102 const std::pair<Vector_t, double>& boundingSphere);
103
105
107 const PART& particle, bool pdead = false);
108
110
111 void calcStat(double Eng);
112 void gatherStatistics();
113
114 void push();
115 void resetTimeStep();
117
118 double T_m; // own time, maybe larger than in the bunch object
119 double dT_m; // dt from bunch
120
121 double mass_m; // mass from bunch (eV)
122 double charge_m; // charge from bunch (elementary charges)
123
124 gsl_rng* rGen_m; // random number generator
125
126 // material parameters
127 std::string material_m; // type of material e.g. aluminum
128 double Z_m; // the atomic number [1]
129 double A_m; // the atomic mass [u]
130 double rho_m; // the volumetric mass density in [g cm^-3]
131 double X0_m; // the radiation length in [m]
132 double I_m; // the mean excitation energy [eV]
133 /*
134 coefficients to fit model to measurement data according to Andersen-Ziegler formulae.
135 see ICRU-49, "Stopping Powers and Ranges for Protons and Alpha Particles",
136 chapter 'Electronic (Collision) Stopping Powers in the Low-Energy Region'
137 */
138 double A1_c;
139 double A2_c;
140 double A3_c;
141 double A4_c;
142 double A5_c;
143 double B1_c;
144 double B2_c;
145 double B3_c;
146 double B4_c;
147 double B5_c;
148
149 // number of particles that enter the material in current step (count for single step)
150 unsigned int bunchToMatStat_m;
151 // number of particles that are stopped by the material in current step
152 unsigned int stoppedPartStat_m;
153 // number of particles that leave the material in current step
154 unsigned int rediffusedStat_m;
155 // total number of particles that are in the material
156 unsigned int totalPartsInMat_m;
157
158 // some statistics
159 double Eavg_m; // average kinetic energy
160 double Emax_m; // maximum kinetic energy
161 double Emin_m; // minimum kinetic energy
162
163 std::vector<PART> locParts_m; // local particles that are in material
164
165 std::unique_ptr<LossDataSink> lossDs_m;
166
169
173};
174
175inline
177 Eavg_m += Eng;
178 if (Emin_m > Eng)
179 Emin_m = Eng;
180 if (Emax_m < Eng)
181 Emax_m = Eng;
182}
183
184inline
186 return T_m;
187}
188
189inline
191 return (element_ref_m->getName() + "_" + name_m);
192}
193
194inline
196 return totalPartsInMat_m;
197}
198
199inline
201 return rediffusedStat_m;
202}
203
204inline
206 return bunchToMatStat_m;
207}
208
209inline
210const std::string ScatteringPhysics::getType() const {
211 return "ScatteringPhysics";
212}
213
214#endif //SCATTERINGPHYSICS_H
const std::string name
virtual const std::string & getName() const
Get element name.
double Qincol
long IDincol
unsigned localID
double DTincol
double Mincol
int Binincol
Vector_t Rincol
Vector_t Pincol
Vector_t Bfincol
Vector_t Efincol
void addBackToBunch(PartBunchBase< double, 3 > *bunch)
virtual const std::string getType() const override
void configureMaterialParameters()
The material of the collimator.
void copyFromBunch(PartBunchBase< double, 3 > *bunch, const std::pair< Vector_t, double > &boundingSphere)
virtual void apply(PartBunchBase< double, 3 > *bunch, const std::pair< Vector_t, double > &boundingSphere) override
void calcStat(double Eng)
virtual size_t getParticlesInMat() override
void addParticleBackToBunch(PartBunchBase< double, 3 > *bunch, const PART &particle, bool pdead=false)
unsigned int bunchToMatStat_m
virtual unsigned getRediffused() override
IpplTimings::TimerRef DegraderLoopTimer_m
void applyRotation(Vector_t &P, Vector_t &R, double xplane, double thetacou)
virtual bool computeEnergyLoss(PartBunchBase< double, 3 > *bunch, Vector_t &P, const double deltat, bool includeFluctuations=true) const override
ScatteringPhysics(const std::string &name, ElementBase *element, std::string &mat, bool enableRutherford, double lowEnergyThr)
void computeInteraction(PartBunchBase< double, 3 > *bunch)
unsigned int totalPartsInMat_m
void applyRandomRotation(Vector_t &P, double theta0)
void computeCoulombScattering(Vector_t &R, Vector_t &P, double dt)
std::vector< PART > locParts_m
virtual unsigned int getNumEntered() override
unsigned int rediffusedStat_m
IpplTimings::TimerRef DegraderApplyTimer_m
unsigned int stoppedPartStat_m
virtual bool stillActive() override
IpplTimings::TimerRef DegraderDestroyTimer_m
virtual std::string getName() override
std::unique_ptr< LossDataSink > lossDs_m
virtual double getTime() override
virtual void print(Inform &os) override
Definition: Inform.h:42
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176