OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
Gaussian.h
Go to the documentation of this file.
1#ifndef IPPL_GAUSSIAN_H
2#define IPPL_GAUSSIAN_H
3
4#include <Kokkos_Random.hpp>
5#include <cmath>
6#include <memory>
7#include "Distribution.h"
8#include "Ippl.h"
9#include "OPALTypes.h"
10#include "SamplingBase.hpp"
11#include "Utilities/Options.h"
12
13using GeneratorPool = typename Kokkos::Random_XorShift64_Pool<>;
14using Dist_t = ippl::random::NormalDistribution<double, 3>;
15
38class Gaussian : public SamplingBase {
39public:
43 IpplTimings::TimerRef samperTimer_m;
44
53 std::shared_ptr<ParticleContainer_t> pc, std::shared_ptr<FieldContainer_t> fc,
54 Distribution_t* opalDist);
55
66 std::shared_ptr<ParticleContainer_t> pc, const Vector_t<double, 3>& sigmaR,
67 const Vector_t<double, 3>& sigmaP, double avrgpz, const Vector_t<double, 3>& cutoffR,
68 bool fix_meanR = true);
75 void generateParticles(size_t& numberOfParticles, Vector_t<double, 3> nr) override;
76
83 void emitParticles(double t, double dt) override;
84
85 void setSigmaR(const Vector_t<double, 3>& sigmaR) { sigmaR_m = sigmaR; }
86
87 void setSigmaP(const Vector_t<double, 3>& sigmaP) { sigmaP_m = sigmaP; }
88
89 void setAvrgpz(double avrgpz) { avrgpz_m = avrgpz; }
90
91 void setCutoffR(const Vector_t<double, 3>& cutoffR) { cutoffR_m = cutoffR; }
92
94 Vector_t<double, 3>& sigmaR, Vector_t<double, 3>& sigmaP, double& avrgpz,
95 Vector_t<double, 3>& cutoffR) const {
96 sigmaR = sigmaR_m;
97 sigmaP = sigmaP_m;
98 avrgpz = avrgpz_m;
99 cutoffR = cutoffR_m;
100 }
101
102 void setFixMeanR(bool fixMeanR) { fixMeanR_m = fixMeanR; }
103
104 void getFixMeanR(bool& fixMeanR) const { fixMeanR = fixMeanR_m; }
105
106private:
110 void initRandomPool();
111
116
122
126 double avrgpz_m;
127
132
136 bool fixMeanR_m = true;
137};
138
139#endif // IPPL_GAUSSIAN_H
ippl::Vector< T, Dim > Vector_t
typename Kokkos::Random_XorShift64_Pool<> GeneratorPool
const int nr
ippl::random::NormalDistribution< double, 3 > Dist_t
Definition FlatTop.cpp:8
typename Kokkos::Random_XorShift64_Pool<> GeneratorPool
Definition Gaussian.h:13
Generating particles following a Gaussian distribution.
Definition Gaussian.h:38
Vector_t< double, 3 > sigmaP_m
Definition Gaussian.h:121
Vector_t< double, 3 > sigmaR_m
Standard deviations for position and momentum distributions.
Definition Gaussian.h:120
void initRandomPool()
Initializes the random number generator pool.
Definition Gaussian.cpp:42
void setCutoffR(const Vector_t< double, 3 > &cutoffR)
Definition Gaussian.h:91
double avrgpz_m
Average momentum in the z-direction.
Definition Gaussian.h:126
void setFixMeanR(bool fixMeanR)
Definition Gaussian.h:102
bool fixMeanR_m
Flag to exactly fix the mean position of particles after sampling.
Definition Gaussian.h:136
void setSigmaR(const Vector_t< double, 3 > &sigmaR)
Definition Gaussian.h:85
void getFixMeanR(bool &fixMeanR) const
Definition Gaussian.h:104
void generateParticles(size_t &numberOfParticles, Vector_t< double, 3 > nr) override
Generates particles with a Gaussian distribution.
Definition Gaussian.cpp:64
Vector_t< double, 3 > cutoffR_m
Cutoff multiplier for position distribution.
Definition Gaussian.h:131
GeneratorPool randPool_m
Pool of random number generators for parallel sampling.
Definition Gaussian.h:115
void getParameters(Vector_t< double, 3 > &sigmaR, Vector_t< double, 3 > &sigmaP, double &avrgpz, Vector_t< double, 3 > &cutoffR) const
Definition Gaussian.h:93
void emitParticles(double t, double dt) override
Time-stepped emission hook for one-shot delayed injection.
Definition Gaussian.cpp:187
void setAvrgpz(double avrgpz)
Definition Gaussian.h:89
void setSigmaP(const Vector_t< double, 3 > &sigmaP)
Definition Gaussian.h:87
IpplTimings::TimerRef samperTimer_m
Timer for performance profiling.
Definition Gaussian.h:43