Code indexing in gitaly is broken and leads to code not being visible to the user. We work on the issue with highest priority.

Skip to content
Snippets Groups Projects

Resolve "Reviewing physics behind particle matter interaction models"

2 files
+ 24
31
Compare changes
  • Side-by-side
  • Inline
Files
2
//
// Class CollimatorPhysics
//
// Defines the collimator physics models
// Defines the collimator physics models
//
// Copyright (c) 2009 - 2020, Bi, Yang, Stachel, Adelmann
// Paul Scherrer Institut, Villigen PSI, Switzerland
@@ -17,7 +16,6 @@
// You should have received a copy of the GNU General Public License
// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
//
#include "Solvers/CollimatorPhysics.hh"
#include "Physics/Physics.h"
#include "Physics/Material.h"
@@ -39,10 +37,10 @@
#include <gsl/gsl_randist.h>
#include <algorithm>
#include <cmath>
#include <iostream>
#include <fstream>
#include <algorithm>
#include <iostream>
#include <sys/time.h>
@@ -152,8 +150,8 @@ CollimatorPhysics::CollimatorPhysics(const std::string& name,
lossDs_m = std::unique_ptr<LossDataSink>(new LossDataSink(getName(), !Options::asciidump));
DegraderApplyTimer_m = IpplTimings::getTimer("DegraderApply");
DegraderLoopTimer_m = IpplTimings::getTimer("DegraderLoop");
DegraderApplyTimer_m = IpplTimings::getTimer("DegraderApply");
DegraderLoopTimer_m = IpplTimings::getTimer("DegraderLoop");
DegraderDestroyTimer_m = IpplTimings::getTimer("DegraderDestroy");
}
@@ -283,15 +281,14 @@ void CollimatorPhysics::computeInteraction() {
}
}
/*
delete absorbed particles
*/
// delete absorbed particles
deleteParticleFromLocalVector();
}
/// Energy Loss: using the Bethe-Bloch equation.
/// Energy straggling: For relatively thick absorbers such that the number of collisions is large,
/// the energy loss distribution is shown to be Gaussian in form.
/// Energy Loss: using the Bethe-Bloch equation.
/// In low-energy region use Anderson-Ziegler fitting (only for protons)
/// Energy straggling: For relatively thick absorbers such that the number of collisions
/// is large, the energy loss distribution is shown to be Gaussian in form.
/// See Particle Physics Booklet, chapter 'Passage of particles through matter' or
/// Review of Particle Physics, DOI: 10.1103/PhysRevD.86.010001, page 329 ff
// -------------------------------------------------------------------------
@@ -345,7 +342,6 @@ bool CollimatorPhysics::computeEnergyLoss(Vector_t& P,
Ekin += deltasrho * dEdx;
}
if (includeFluctuations) {
double sigma_E = std::sqrt(K * massElectron_keV * rho_m * (Z_m / A_m) * deltas * m2cm);
Ekin += gsl_ran_gaussian(rGen_m, sigma_E);
@@ -361,9 +357,9 @@ bool CollimatorPhysics::computeEnergyLoss(Vector_t& P,
// Implement the rotation in 2 dimensions here
// For details see:
// J. Beringer et al. (Particle Data Group), "Passage of particles through matter"
// Phys. Rev. D 86, 010001 (2012),
// For details see: J. Beringer et al. (Particle Data Group),
// "Passage of particles through matter", Phys. Rev. D 86, 010001 (2012)
// -------------------------------------------------------------------------
void CollimatorPhysics::applyRotation(Vector_t& P,
Vector_t& R,
double shift,
@@ -488,9 +484,7 @@ void CollimatorPhysics::addBackToBunch(PartBunchBase<double, 3>* bunch) {
}
}
/*
delete particles that went to the bunch
*/
// delete particles that went to the bunch
deleteParticleFromLocalVector();
}
@@ -596,8 +590,8 @@ namespace {
void CollimatorPhysics::deleteParticleFromLocalVector() {
/*
the particle to be deleted (label < 0) are all at the end of
the vector.
the particle to be deleted (label < 0) are all
at the end of the vector.
*/
sort(locParts_m.begin(), locParts_m.end(), myCompF);
Loading