OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
FlatTop.h
Go to the documentation of this file.
1
6#ifndef IPPL_FLAT_TOP_H
7#define IPPL_FLAT_TOP_H
8
9#include <Kokkos_Random.hpp>
10#include <cmath>
11#include <memory>
12#include "Distribution.h"
13#include "Ippl.h"
14#include "OPALTypes.h"
15#include "SamplingBase.hpp"
16#include "Utilities/Options.h"
17
18using GeneratorPool = typename Kokkos::Random_XorShift64_Pool<>;
19using Dist_t = ippl::random::NormalDistribution<double, 3>;
20
33class FlatTop : public SamplingBase {
34public:
41 FlatTop(std::shared_ptr<ParticleContainer_t> pc, std::shared_ptr<FieldContainer_t> fc,
42 Distribution_t* opalDist);
43
56 FlatTop(std::shared_ptr<ParticleContainer_t> pc, std::shared_ptr<FieldContainer_t> fc,
57 bool emitting, double sigmaTFall, double sigmaTRise, Vector_t<double, 3> cutoff,
58 double tPulseLengthFWHM, Vector_t<double, 3> sigmaR);
59
70 FlatTop(std::shared_ptr<ParticleContainer_t> pc, bool emitting, double sigmaTFall,
71 double sigmaTRise, Vector_t<double, 3> cutoff, double tPulseLengthFWHM,
72 Vector_t<double, 3> sigmaR);
73
79 void testNumEmitParticles(size_type nsteps, double dt) override;
80
86 void testEmitParticles(size_type nsteps, double dt) override;
87
92 void setWithDomainDecomp(bool withDomainDecomp) override;
93
97 void setDistArea(double distArea) { distArea_m = distArea; }
98
103 double getDistArea() { return distArea_m; }
104
108 double getEmissionTime() { return emissionTime_m; }
109
111 bool isEmissionDone(double t) const override { return (t - t0_m) >= emissionTime_m; }
112
113private:
114 using size_type = ippl::detail::size_type;
118 double distArea_m;
122 double fallTime_m;
123 double riseTime_m;
131
136 static size_t determineRandInit();
137
142 void setParameters(Distribution_t* opalDist);
143
144public:
157 void generateUniformDisk(size_type nlocal, size_t nNew, double dt);
158
164
170 void generateParticles(size_t& numberOfParticles, Vector_t<double, 3> nr) override;
171
177 double FlatTopProfile(double t);
178
184 size_t computeNlocalUniformly(size_t nglobal);
185
194 double integrateTrapezoidal(double x1, double x2, double y1, double y2);
195
200 void initDomainDecomp(double BoxIncr) override;
201
208 size_type countEnteringParticlesPerRank(double t0, double tf);
209
214 void allocateParticles(size_t numberOfParticles);
215
221 void emitParticles(double t, double dt) override;
222
224 bool emitting, double sigmaTFall, double sigmaTRise, Vector_t<double, 3> cutoff,
225 double tPulseLengthFWHM, Vector_t<double, 3> sigmaR);
226};
227
228#endif // IPPL_FLAT_TOP_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 FlatTop.h:18
Implements the sampling method for the flat-top distribution. x and y coordinates are uniformly distr...
Definition FlatTop.h:33
void setNr(Vector_t< double, 3 > nr)
Sets the number of grid points per direction.
Definition FlatTop.cpp:181
void setDistArea(double distArea)
Sets the total area of the distribution.
Definition FlatTop.h:97
void generateUniformDisk(size_type nlocal, size_t nNew, double dt)
Generates particles (x,y) uniformly on a disk distribution.
Definition FlatTop.cpp:88
void testEmitParticles(size_type nsteps, double dt) override
Tests particle emission over a given number of steps.
Definition FlatTop.cpp:440
double distArea_m
Total area of the flattop distribution.
Definition FlatTop.h:118
Vector_t< double, 3 > sigmaR_m
Definition FlatTop.h:130
double FlatTopProfile(double t)
Computes the flat-top profile value at a given time.
Definition FlatTop.cpp:204
bool emitting_m
Flag for particle emission status.
Definition FlatTop.h:124
Vector_t< double, 3 > nr_m
Number of grid points per direction.
Definition FlatTop.h:128
static size_t determineRandInit()
Determines the random seed initialization.
Definition FlatTop.cpp:35
void allocateParticles(size_t numberOfParticles)
Allocates memory for a given number of particles.
Definition FlatTop.cpp:350
double emissionTime_m
Total emission time.
Definition FlatTop.h:127
double sigmaTRise_m
Standard deviation for rise time profile.
Definition FlatTop.h:120
size_type totalN_m
Total number of particles.
Definition FlatTop.h:125
GeneratorPool rand_pool_m
Random number generator pool.
Definition FlatTop.h:115
double sigmaTFall_m
Standard deviation for fall time profile.
Definition FlatTop.h:119
bool withDomainDecomp_m
Flag for domain decomposition.
Definition FlatTop.h:126
void initDomainDecomp(double BoxIncr) override
Initializes the domain decomposition.
Definition FlatTop.cpp:242
Vector_t< double, 3 > hr_m
Grid spacing.
Definition FlatTop.h:129
double integrateTrapezoidal(double x1, double x2, double y1, double y2)
Integrates using the trapezoidal rule.
Definition FlatTop.cpp:238
double fallTime_m
Time duration for the fall phase.
Definition FlatTop.h:122
size_t computeNlocalUniformly(size_t nglobal)
Computes the local number of particles uniformly distributed among ranks.
Definition FlatTop.cpp:222
double normalizedFlankArea_m
Normalized area of the distribution flanks.
Definition FlatTop.h:117
void setWithDomainDecomp(bool withDomainDecomp) override
Sets whether to use domain decomposition.
Definition FlatTop.cpp:33
double getEmissionTime()
Returns the total emission time.
Definition FlatTop.h:108
void emitParticles(double t, double dt) override
Emits new particles within a given time interval.
Definition FlatTop.cpp:360
void setInternalVariables(bool emitting, double sigmaTFall, double sigmaTRise, Vector_t< double, 3 > cutoff, double tPulseLengthFWHM, Vector_t< double, 3 > sigmaR)
Definition FlatTop.cpp:58
double riseTime_m
Time duration for the rise phase.
Definition FlatTop.h:123
void generateParticles(size_t &numberOfParticles, Vector_t< double, 3 > nr) override
Generates particles with a given number and grid configuration.
Definition FlatTop.cpp:183
void setParameters(Distribution_t *opalDist)
Sets distribution parameters.
Definition FlatTop.cpp:47
void testNumEmitParticles(size_type nsteps, double dt) override
Tests the number of emitted particles over a given number of steps.
Definition FlatTop.cpp:394
Vector_t< double, 3 > cutoffR_m
Cutoff radius.
Definition FlatTop.h:121
double flattopTime_m
Time duration of when the time profile is flat.
Definition FlatTop.h:116
size_type countEnteringParticlesPerRank(double t0, double tf)
Counts the number of particles entering per rank in a given time interval.
Definition FlatTop.cpp:263
ippl::detail::size_type size_type
Definition FlatTop.h:114
bool isEmissionDone(double t) const override
Whether this sampler has finished all emission (no more particles will be created).
Definition FlatTop.h:111
double getDistArea()
Returns the total area of the distribution.
Definition FlatTop.h:103