OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
ScatteringPhysics.cpp
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//
21
26#include "AbsBeamline/Drift.h"
27#include "AbsBeamline/SBend.h"
28#include "AbsBeamline/RBend.h"
30#include "Physics/Material.h"
32#include "Physics/Physics.h"
33#include "Physics/Units.h"
35#include "Utilities/Options.h"
37#include "Utilities/Util.h"
38#include "Utilities/Timer.h"
39
40#include "Utility/Inform.h"
41
42#include <gsl/gsl_randist.h>
43
44#include <algorithm>
45#include <cmath>
46#include <filesystem>
47#include <fstream>
48#include <iostream>
49
50#include <sys/time.h>
51
52namespace {
53 struct DegraderInsideTester: public InsideTester {
54 explicit DegraderInsideTester(ElementBase* el) {
55 deg_m = static_cast<Degrader*>(el);
56 }
57 bool checkHit(const Vector_t& R) override {
58 return deg_m->isInside(R);
59 }
60 private:
61 Degrader* deg_m;
62 };
63
64 struct CollimatorInsideTester: public InsideTester {
65 explicit CollimatorInsideTester(ElementBase* el) {
66 col_m = static_cast<CCollimator*>(el);
67 }
68 bool checkHit(const Vector_t& R) override {
69 return col_m->checkPoint(R(0), R(1));
70 }
71 private:
72 CCollimator* col_m;
73 };
74
75 struct FlexCollimatorInsideTester: public InsideTester {
76 explicit FlexCollimatorInsideTester(ElementBase* el) {
77 col_m = static_cast<FlexibleCollimator*>(el);
78 }
79 bool checkHit(const Vector_t& R) override {
80 return col_m->isStopped(R);
81 }
82 private:
83 FlexibleCollimator* col_m;
84 };
85
86 constexpr long double operator"" _keV(long double value) { return value; }
87 constexpr long double operator"" _MeV(long double value) { return value * 1e3; }
88 constexpr long double operator"" _GeV(long double value) { return value * 1e6; }
89}
90
92 ElementBase* element,
93 std::string& material,
94 bool enableRutherford,
95 double lowEnergyThr):
97 T_m(0.0),
98 dT_m(0.0),
99 mass_m(0.0),
100 charge_m(0.0),
101 material_m(material),
102 Z_m(0),
103 A_m(0.0),
104 rho_m(0.0),
105 X0_m(0.0),
106 I_m(0.0),
107 A1_c(0.0),
108 A2_c(0.0),
109 A3_c(0.0),
110 A4_c(0.0),
111 A5_c(0.0),
112 B1_c(0.0),
113 B2_c(0.0),
114 B3_c(0.0),
115 B4_c(0.0),
116 B5_c(0.0),
121 Eavg_m(0.0),
122 Emax_m(0.0),
123 Emin_m(0.0),
124 enableRutherford_m(enableRutherford),
125 lowEnergyThr_m(lowEnergyThr)
126{
127
128 gsl_rng_env_setup();
129 rGen_m = gsl_rng_alloc(gsl_rng_default);
130
131 size_t mySeed = Options::seed;
132
133 if (Options::seed == -1) {
134 struct timeval tv;
135 gettimeofday(&tv,0);
136 mySeed = tv.tv_sec + tv.tv_usec;
137 }
138
139 gsl_rng_set(rGen_m, mySeed + Ippl::myNode());
140
142
143 ElementType collshape = element_ref_m->getType();
144 switch (collshape) {
146 hitTester_m.reset(new DegraderInsideTester(element_ref_m));
147 break;
149 hitTester_m.reset(new CollimatorInsideTester(element_ref_m));
150 break;
152 hitTester_m.reset(new FlexCollimatorInsideTester(element_ref_m));
153 break;
154 default:
155 throw GeneralClassicException("ScatteringPhysics::ScatteringPhysics",
156 "Unsupported element type");
157 }
158
159 lossDs_m = std::unique_ptr<LossDataSink>(new LossDataSink(getName(), !Options::asciidump));
160
164}
165
167 locParts_m.clear();
168 lossDs_m->save();
169 if (rGen_m)
170 gsl_rng_free(rGen_m);
171}
172
173
175// ------------------------------------------------------------------------
178 Z_m = material->getAtomicNumber();
179 A_m = material->getAtomicMass();
180 rho_m = material->getMassDensity();
181 X0_m = material->getRadiationLength();
182 I_m = material->getMeanExcitationEnergy();
183 A1_c = material->getStoppingPowerFitCoefficients(Physics::Material::A1);
184 A2_c = material->getStoppingPowerFitCoefficients(Physics::Material::A2);
185 A3_c = material->getStoppingPowerFitCoefficients(Physics::Material::A3);
186 A4_c = material->getStoppingPowerFitCoefficients(Physics::Material::A4);
187 A5_c = material->getStoppingPowerFitCoefficients(Physics::Material::A5);
188 B1_c = material->getStoppingPowerFitCoefficients(Physics::Material::B1);
189 B2_c = material->getStoppingPowerFitCoefficients(Physics::Material::B2);
190 B3_c = material->getStoppingPowerFitCoefficients(Physics::Material::B3);
191 B4_c = material->getStoppingPowerFitCoefficients(Physics::Material::B4);
192 B5_c = material->getStoppingPowerFitCoefficients(Physics::Material::B5);
193}
194
196 const std::pair<Vector_t, double>& boundingSphere) {
198
199 /*
200 Particles that have entered material are flagged as Bin[i] == -1.
201
202 Flagged particles are copied to a local structure within Scattering Physics locParts_m.
203
204 Particles in that structure will be pushed in the material and either come
205 back to the bunch or will be fully stopped in the material.
206 */
207
208 ParticleType pType = bunch->getPType();
209 if (pType != ParticleType::PROTON &&
210 pType != ParticleType::DEUTERON &&
211 pType != ParticleType::HMINUS &&
212 pType != ParticleType::MUON &&
213 pType != ParticleType::H2P &&
214 pType != ParticleType::ALPHA) {
215
217 "ScatteringPhysics::apply",
218 "Particle " + ParticleProperties::getParticleTypeString(pType) +
219 " is not supported for scattering interactions!");
220 }
221
222 Eavg_m = 0.0;
223 Emax_m = 0.0;
224 Emin_m = 0.0;
225
230
231 dT_m = bunch->getdT();
232 T_m = bunch->getT();
233 mass_m = bunch->getM();
234 charge_m = bunch->getQ();
235
236 bool onlyOneLoopOverParticles = !(allParticleInMat_m);
237
238 do {
240 push();
242
243 // if we are not looping copy newly arrived particles
244 if (onlyOneLoopOverParticles) {
245 copyFromBunch(bunch, boundingSphere);
246 }
247 addBackToBunch(bunch);
248
249 computeInteraction(bunch);
250
251 push();
253
255
256 T_m += dT_m; // update local time
257
259
260 if (allParticleInMat_m) {
261 onlyOneLoopOverParticles = (rediffusedStat_m > 0 || totalPartsInMat_m <= 0);
262 } else {
263 onlyOneLoopOverParticles = true;
264 }
265 } while (onlyOneLoopOverParticles == false);
266
268}
269
271 /*
272 Do physics if
273 -- correct type of particle
274 -- particle not stopped (locParts_m[i].label != -1.0)
275
276 Absorbed particle i: locParts_m[i].label = -1.0;
277 */
278 for (size_t i = 0; i < locParts_m.size(); ++i) {
279 if (locParts_m[i].label != -1) {
280 Vector_t& R = locParts_m[i].Rincol;
281 Vector_t& P = locParts_m[i].Pincol;
282 double& dt = locParts_m[i].DTincol;
283
284 if (hitTester_m->checkHit(R)) {
285 bool pdead = computeEnergyLoss(bunch, P, dt);
286 if (!pdead) {
287 /*
288 Now scatter and transport particle in material.
289 The checkHit call just above will detect if the
290 particle is rediffused from the material into vacuum.
291 */
292
293 computeCoulombScattering(R, P, dt);
294 } else {
295 // The particle is stopped in the material, set label to -1
296 locParts_m[i].label = -1.0;
298 lossDs_m->addParticle(OpalParticle(locParts_m[i].IDincol,
299 R, P, T_m,
300 locParts_m[i].Qincol, locParts_m[i].Mincol));
301
302 if (OpalData::getInstance()->isInOPALCyclMode()) {
303 // OpalCycl performs particle deletion in ParallelCyclotronTracker.
304 // Particles lost by scattering have to be returned to the bunch
305 // with a negative Bin attribute to avoid miscounting of particles.
306 // Only the minimal number of attributes are fixed because the
307 // particle is marked for deletion (Bin<0)
308
309 addParticleBackToBunch(bunch, locParts_m[i], pdead);
310 }
311 }
312 }
313 }
314 }
315
316 // delete absorbed particles
318}
319
326// -------------------------------------------------------------------------
328 Vector_t& P,
329 const double deltat,
330 bool includeFluctuations) const {
331
332 ParticleType pType = bunch->getPType();
333
334 const double mass_keV = mass_m * Units::eV2keV;
335 constexpr double massElectron_keV = Physics::m_e * Units::GeV2keV;
336
337 constexpr double K = (4.0 * Physics::pi * Physics::Avo * Physics::r_e
338 * Units::m2cm * Physics::r_e * Units::m2cm * massElectron_keV);
339
340 double gamma = Util::getGamma(P);
341 const double gammaSqr = std::pow(gamma, 2);
342 const double betaSqr = 1.0 - 1.0 / gammaSqr;
343 double beta = std::sqrt(betaSqr);
344 double Ekin_keV = (gamma - 1) * mass_keV;
345 double dEdx = 0.0;
346 double epsilon = 0.0;
347
348 const double deltas = deltat * beta * Physics::c;
349 const double deltasrho = deltas * Units::m2cm * rho_m;
350
351 const double massRatio = massElectron_keV / mass_keV;
352 double Tmax = (2.0 * massElectron_keV * betaSqr * gammaSqr /
353 (std::pow(gamma + massRatio, 2) - (gammaSqr - 1.0)));
354
355 if (pType != ParticleType::ALPHA) {
356
357 if (Ekin_keV >= 0.6_MeV && Ekin_keV < 10.0_GeV) {
358 dEdx = (-K * std::pow(charge_m, 2) * Z_m / (A_m * betaSqr) *
359 (0.5 * std::log(2 * massElectron_keV * betaSqr * gammaSqr * Tmax / std::pow(I_m * Units::eV2keV, 2)) - betaSqr));
360 } else if (pType == ParticleType::PROTON && Ekin_keV < 0.6_MeV) {
361 constexpr double massProton_amu = Physics::m_p / Physics::amu;
362 const double Ts = Ekin_keV / massProton_amu;
363 if (Ekin_keV > 10.0_keV) {
364 const double epsilon_low = A2_c * std::pow(Ts, 0.45);
365 const double epsilon_high = (A3_c / Ts) * std::log(1 + (A4_c / Ts) + (A5_c * Ts));
366 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
367 } else if (Ekin_keV > 1.0_keV) {
368 epsilon = A1_c * std::pow(Ts, 0.5);
369 }
370 dEdx = -epsilon / (1e18 * (A_m / Physics::Avo));
371 } else {
372 INFOMSG(level4 << "Particle energy out of the valid range "
373 "for energy loss calculation." << endl);
374 }
375 } else {
376 if (Ekin_keV > 10.0_MeV && Ekin_keV < 1.0_GeV) {
377 dEdx = (-K * std::pow(charge_m, 2) * Z_m / (A_m * betaSqr) *
378 (0.5 * std::log(2 * massElectron_keV * betaSqr * gammaSqr * Tmax / std::pow(I_m * Units::eV2keV, 2)) - betaSqr));
379 } else if (Ekin_keV > 1.0_keV && Ekin_keV <= 10.0_MeV) {
380 const double T = Ekin_keV * Units::keV2MeV;
381 const double epsilon_low = B1_c * std::pow(T * Units::MeV2keV, B2_c);
382 const double epsilon_high = (B3_c / T) * std::log(1 + (B4_c / T) + (B5_c * T));
383 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
384 dEdx = -epsilon / (1e18 * (A_m / Physics::Avo));
385 } else {
386 INFOMSG(level4 << "Particle energy out of the valid range "
387 "for energy loss calculation." << endl);
388 }
389 }
390
391 Ekin_keV += deltasrho * dEdx;
392
393 if (includeFluctuations) {
394 double sigma_E = std::sqrt(K * massElectron_keV * rho_m * (Z_m / A_m) * deltas * Units::m2cm);
395 Ekin_keV += gsl_ran_gaussian(rGen_m, sigma_E);
396 }
397
398 gamma = Ekin_keV / mass_keV + 1.0;
399 beta = std::sqrt(1.0 - 1.0 / std::pow(gamma, 2));
400 P = gamma * beta * P / euclidean_norm(P);
401
402 bool stopped = (Ekin_keV < lowEnergyThr_m * Units::MeV2keV || dEdx > 0);
403 return stopped;
404}
405
406
407// Implement the rotation in 2 dimensions here
408// For details see: J. Beringer et al. (Particle Data Group),
409// "Passage of particles through matter", Phys. Rev. D 86, 010001 (2012)
410// -------------------------------------------------------------------------
412 Vector_t& R,
413 double shift,
414 double thetacou) {
415 // Calculate the angle between the transverse and longitudinal component of the momentum
416 double Psixz = std::fmod(std::atan2(P(0), P(2)) + Physics::two_pi, Physics::two_pi);
417
418 R(0) = R(0) + shift * std::cos(Psixz);
419 R(2) = R(2) - shift * std::sin(Psixz);
420
421 // Apply the rotation about the random angle thetacou
422 double Px = P(0);
423 P(0) = Px * std::cos(thetacou) + P(2) * std::sin(thetacou);
424 P(2) = -Px * std::sin(thetacou) + P(2) * std::cos(thetacou);
425}
426
428
429 double thetaru = 2.5 / std::sqrt(gsl_rng_uniform(rGen_m)) * 2.0 * theta0;
430 double phiru = Physics::two_pi * gsl_rng_uniform(rGen_m);
431
432 double normPtrans = std::sqrt(P(0) * P(0) + P(1) * P(1));
433 double Theta = std::atan(normPtrans / std::abs(P(2)));
434 double CosT = std::cos(Theta);
435 double SinT = std::sin(Theta);
436
437 Vector_t X(std::cos(phiru)*std::sin(thetaru),
438 std::sin(phiru)*std::sin(thetaru),
439 std::cos(thetaru));
440 X *= euclidean_norm(P);
441
442 Vector_t W(-P(1), P(0), 0.0);
443 W = W / normPtrans;
444
445 // Rodrigues' formula for a rotation about W by angle Theta
446 P = X * CosT + cross(W, X) * SinT + W * dot(W, X) * (1.0 - CosT);
447}
448
451//--------------------------------------------------------------------------
453 Vector_t& P,
454 double dt) {
455
456 constexpr double sqrtThreeInv = 0.57735026918962576451; // sqrt(1.0 / 3.0)
457 const double normP = euclidean_norm(P);
458 const double beta = euclidean_norm(Util::getBeta(P));
459 const double deltas = dt * beta * Physics::c;
460 const double theta0 = (13.6e6 / (beta * normP * mass_m) *
461 charge_m * std::sqrt(deltas / X0_m) *
462 (1.0 + 0.038 * std::log(deltas / X0_m)));
463
464 double phi = Physics::two_pi * gsl_rng_uniform(rGen_m);
465 for (unsigned int i = 0; i < 2; ++ i) {
466 CoordinateSystemTrafo randomTrafo(R, Quaternion(cos(phi), 0, 0, sin(phi)));
467 P = randomTrafo.rotateTo(P);
468 R = Vector_t(0.0); // corresponds to randomTrafo.transformTo(R);
469
470 double z1 = gsl_ran_gaussian(rGen_m, 1.0);
471 double z2 = gsl_ran_gaussian(rGen_m, 1.0);
472
473 while(std::abs(z2) > 3.5) {
474 z1 = gsl_ran_gaussian(rGen_m, 1.0);
475 z2 = gsl_ran_gaussian(rGen_m, 1.0);
476 }
477
478 double thetacou = z2 * theta0;
479 double shift = (z1 * sqrtThreeInv + z2) * deltas * theta0 * 0.5;
480 applyRotation(P, R, shift, thetacou);
481
482 P = randomTrafo.rotateFrom(P);
483 R = randomTrafo.transformFrom(R);
484
485 phi += 0.5 * Physics::pi;
486 }
487
488 if (enableRutherford_m && gsl_rng_uniform(rGen_m) < 0.0047) {
489 applyRandomRotation(P, theta0);
490 }
491}
492
494
495 const size_t nL = locParts_m.size();
496 if (nL == 0) return;
497
498 const double elementLength = element_ref_m->getElementLength();
499
500 for (size_t i = 0; i < nL; ++ i) {
501 Vector_t& R = locParts_m[i].Rincol;
502
503 if ( (OpalData::getInstance()->isInOPALTMode() && R[2] >= elementLength) ||
504 (OpalData::getInstance()->isInOPALCyclMode() && !(hitTester_m->checkHit(R))) ) {
505
507
508 /*
509 This particle is back to the bunch, by setting the
510 label to -1.0 the particle will be deleted.
511 */
512 locParts_m[i].label = -1.0;
513
515 }
516 }
517
518 // delete particles that went to the bunch
520}
521
523 const PART& particle, bool pdead) {
524
525 unsigned int numLocalParticles = bunch->getLocalNum();
526
527 bunch->createWithID(particle.IDincol);
528
529 if (!pdead) {
530 bunch->Bin[numLocalParticles] = 1;
531 } else {
532 bunch->Bin[numLocalParticles] = -1;
533 }
534
535 bunch->R[numLocalParticles] = particle.Rincol;
536 bunch->P[numLocalParticles] = particle.Pincol;
537 bunch->Q[numLocalParticles] = particle.Qincol;
538 bunch->M[numLocalParticles] = particle.Mincol;
539 bunch->Bf[numLocalParticles] = 0.0;
540 bunch->Ef[numLocalParticles] = 0.0;
541 bunch->dt[numLocalParticles] = dT_m;
542}
543
544
546 const std::pair<Vector_t, double>& boundingSphere)
547{
548 const size_t nL = bunch->getLocalNum();
549 if (nL == 0) return;
550
552 double zmax = boundingSphere.first(2) + boundingSphere.second;
553 double zmin = boundingSphere.first(2) - boundingSphere.second;
554 if (zmax < 0.0 || zmin > element_ref_m->getElementLength()) {
556 return;
557 }
558
559 size_t ne = 0;
560 std::set<size_t> partsToDel;
561 for (size_t i = 0; i < nL; ++ i) {
562 if ((bunch->Bin[i] == -1 || bunch->Bin[i] == 1) &&
563 hitTester_m->checkHit(bunch->R[i]))
564 {
565 double tau = 1.0;
566 if (OpalData::getInstance()->isInOPALTMode()) {
567 // The z-coordinate is only Opal-T mode the longitudinal coordinate and
568 // the case when elements with ScatteringPhysics solver are closer than one
569 // time step needs to be handled, which isn't done yet in Opal-cycl.
570
571 // Adjust the time step for those particles that enter the material
572 // such that it corresponds to the time needed to reach the current
573 // location form the edge of the material. Only use this time step
574 // for the computation of the interaction with the material, not for
575 // the integration of the particles. This will ensure that the momenta
576 // of all particles are reduced by approximately the same amount in
577 // computeEnergyLoss.
578
579 double betaz = bunch->P[i](2) / Util::getGamma(bunch->P[i]);
580 double stepWidth = betaz * Physics::c * bunch->dt[i];
581 tau = bunch->R[i](2) / stepWidth;
582
583 PAssert_LE(tau, 1.0);
584 PAssert_GE(tau, 0.0);
585 }
586
587 PART x;
588 x.localID = i;
589 x.DTincol = bunch->dt[i] * tau;
590 x.IDincol = bunch->ID[i];
591 x.Binincol = bunch->Bin[i];
592 x.Rincol = bunch->R[i];
593 x.Pincol = bunch->P[i];
594 x.Qincol = bunch->Q[i];
595 x.Mincol = bunch->M[i];
596 x.Bfincol = bunch->Bf[i];
597 x.Efincol = bunch->Ef[i];
598 x.label = 0; // alive in matter
599
600 locParts_m.push_back(x);
601 ++ne;
603
604 partsToDel.insert(i);
605 }
606 }
607
608 for (auto it = partsToDel.begin(); it != partsToDel.end(); ++ it) {
609 bunch->destroy(1, *it);
610 }
611
612 if (ne > 0) {
613 bunch->performDestroy(true);
614 }
615
617}
618
620 Inform::FmtFlags_t ff = msg.flags();
621 if (totalPartsInMat_m > 0 ||
622 bunchToMatStat_m > 0 ||
623 rediffusedStat_m > 0 ||
624 stoppedPartStat_m > 0) {
625
626 OPALTimer::Timer time;
627 msg << level2
628 << "--- ScatteringPhysics ---\n"
629 << "Name: " << name_m << " - "
630 << "Material: " << material_m << " - "
631 << "Element: " << element_ref_m->getName() << "\n"
632 << "Particle Statistics @ " << time.time() << "\n"
633 << std::setw(21) << "entered: " << Util::toStringWithThousandSep(bunchToMatStat_m) << "\n"
634 << std::setw(21) << "rediffused: " << Util::toStringWithThousandSep(rediffusedStat_m) << "\n"
635 << std::setw(21) << "stopped: " << Util::toStringWithThousandSep(stoppedPartStat_m) << "\n"
636 << std::setw(21) << "total in material: " << Util::toStringWithThousandSep(totalPartsInMat_m)
637 << endl;
638 }
639 msg.flags(ff);
640}
641
645
646namespace {
647 bool myCompF(PART x, PART y) {
648 return x.label > y.label;
649 }
650}
651
653 /*
654 the particle to be deleted (label < 0) are all
655 at the end of the vector.
656 */
657 sort(locParts_m.begin(), locParts_m.end(), myCompF);
658
659 // find start of particles to delete
660 std::vector<PART>::iterator inv = locParts_m.begin();
661
662 for (; inv != locParts_m.end(); ++inv) {
663 if ((*inv).label == -1)
664 break;
665 }
666 locParts_m.erase(inv, locParts_m.end());
667 locParts_m.resize(inv - locParts_m.begin());
668
669 // update statistics
670 if (!locParts_m.empty()) {
671 Eavg_m /= locParts_m.size();
672 Emin_m /= locParts_m.size();
673 Emax_m /= locParts_m.size();
674 }
675}
676
678 for (size_t i = 0; i < locParts_m.size(); ++i) {
679 Vector_t& R = locParts_m[i].Rincol;
680 Vector_t& P = locParts_m[i].Pincol;
681 double gamma = Util::getGamma(P);
682
683 R += 0.5 * dT_m * Physics::c * P / gamma;
684 }
685}
686
687// adjust the time step for those particles that will rediffuse to
688// vacuum such that it corresponds to the time needed to reach the
689// end of the material. Only use this time step for the computation
690// of the interaction with the material, not for the integration of
691// the particles. This will ensure that the momenta of all particles
692// are reduced by approximately the same amount in computeEnergyLoss.
694 const double elementLength = element_ref_m->getElementLength();
695
696 for (size_t i = 0; i < locParts_m.size(); ++i) {
697 Vector_t& R = locParts_m[i].Rincol;
698 Vector_t& P = locParts_m[i].Pincol;
699 double& dt = locParts_m[i].DTincol;
700 double gamma = Util::getGamma(P);
701 Vector_t stepLength = dT_m * Physics::c * P / gamma;
702
703 if ( R(2) < elementLength &&
704 R(2) + stepLength(2) > elementLength &&
705 OpalData::getInstance()->isInOPALTMode() ) {
706
707 // particle is likely to leave material
708 double distance = elementLength - R(2);
709 double tau = distance / stepLength(2);
710
711 PAssert_LE(tau, 1.0);
712 PAssert_GE(tau, 0.0);
713
714 dt *= tau;
715 }
716 }
717}
718
720 for (size_t i = 0; i < locParts_m.size(); ++i) {
721 double& dt = locParts_m[i].DTincol;
722 dt = dT_m;
723 }
724}
725
727
728 unsigned int locPartsInMat;
729 locPartsInMat = locParts_m.size();
730
731 constexpr unsigned short numItems = 4;
732 unsigned int partStatistics[numItems] = {locPartsInMat,
736
737 allreduce(&(partStatistics[0]), numItems, std::plus<unsigned int>());
738
739 totalPartsInMat_m = partStatistics[0];
740 bunchToMatStat_m = partStatistics[1];
741 rediffusedStat_m = partStatistics[2];
742 stoppedPartStat_m = partStatistics[3];
743}
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
Definition Vector3D.cpp:111
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition Vector3D.cpp:118
ElementType
Definition ElementBase.h:88
struct @313225227106163032040204343033076100017014020135 PART
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition TpsMath.h:129
Tps< T > sin(const Tps< T > &x)
Sine.
Definition TpsMath.h:111
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition Vector.h:243
void allreduce(const T *input, T *output, int count, Op op)
PETE_TBTree< OpNE, Index::PETE_Expr_t, PETE_Scalar< double > > ne(const Index &idx, double x)
#define X(arg)
Definition fftpack.cpp:106
Inform & level2(Inform &inf)
Definition Inform.cpp:46
Inform & level4(Inform &inf)
Definition Inform.cpp:48
Inform & endl(Inform &inf)
Definition Inform.cpp:42
#define PAssert_LE(a, b)
Definition PAssert.h:107
#define PAssert_GE(a, b)
Definition PAssert.h:109
#define INFOMSG(msg)
Definition IpplInfo.h:348
const std::string name
constexpr double amu
The atomic mass unit energy equivalent in GeV.
Definition Physics.h:75
constexpr double two_pi
The value of.
Definition Physics.h:33
constexpr double Avo
The Avogadro's number.
Definition Physics.h:57
constexpr double m_p
The proton rest mass in GeV.
Definition Physics.h:90
constexpr double m_e
The electron rest mass in GeV.
Definition Physics.h:78
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
constexpr double pi
The value of.
Definition Physics.h:30
constexpr double r_e
The classical electron radius in m.
Definition Physics.h:81
constexpr double GeV2keV
Definition Units.h:92
constexpr double eV2keV
Definition Units.h:89
constexpr double m2cm
Definition Units.h:32
constexpr double keV2MeV
Definition Units.h:101
constexpr double MeV2keV
Definition Units.h:98
bool asciidump
Definition Options.cpp:85
int seed
The current random seed.
Definition Options.cpp:37
Vector_t getBeta(Vector_t p)
Definition Util.h:51
double getGamma(Vector_t p)
Definition Util.h:46
std::string toStringWithThousandSep(T value, char sep='\'')
Definition Util.h:254
ParticleAttrib< Vector_t > Ef
ParticlePos_t & R
ParticleAttrib< int > Bin
double getQ() const
Access to reference data.
size_t getLocalNum() const
ParticleAttrib< double > M
ParticleAttrib< Vector_t > P
double getdT() const
ParticleAttrib< double > Q
ParticleType getPType() const
void createWithID(unsigned id)
ParticleAttrib< double > dt
void performDestroy(bool updateLocalNum=false)
void destroy(size_t M, size_t I, bool doNow=false)
ParticleAttrib< Vector_t > Bf
ParticleIndex_t & ID
double getM() const
double getT() const
static OpalData * getInstance()
Definition OpalData.cpp:196
Vector_t rotateFrom(const Vector_t &r) const
Vector_t transformFrom(const Vector_t &r) const
Vector_t rotateTo(const Vector_t &r) const
static std::shared_ptr< Material > getMaterial(const std::string &name)
Definition Material.cpp:53
static std::string getParticleTypeString(const ParticleType &type)
ParticleMatterInteractionHandler(std::string name, ElementBase *elref)
bool allParticleInMat_m
if all particles are in matter stay inside the particle matter interaction
void addBackToBunch(PartBunchBase< double, 3 > *bunch)
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 addParticleBackToBunch(PartBunchBase< double, 3 > *bunch, const PART &particle, bool pdead=false)
unsigned int bunchToMatStat_m
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
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 void print(Inform &os) override
std::string time() const
Return time.
std::ios_base::fmtflags FmtFlags_t
Definition Inform.h:99
FmtFlags_t flags() const
Definition Inform.h:106
static int myNode()
Definition IpplInfo.cpp:691
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)
Vektor< double, 3 > Vector_t