OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
FlatTop.h
Go to the documentation of this file.
1
5
6#ifndef IPPL_FLAT_TOP_H
7#define IPPL_FLAT_TOP_H
8
9#include "Distribution.h"
10#include "SamplingBase.hpp"
11#include <Kokkos_Random.hpp>
12#include "Ippl.h"
13#include "Utilities/Options.h"
14#include "OPALTypes.h"
15#include <memory>
16#include <cmath>
17
21using GeneratorPool = typename Kokkos::Random_XorShift64_Pool<>;
22using Dist_t = ippl::random::NormalDistribution<double, 3>;
23
36class FlatTop : public SamplingBase {
37public:
44 FlatTop(std::shared_ptr<ParticleContainer_t> &pc, std::shared_ptr<FieldContainer_t> &fc, std::shared_ptr<Distribution_t> &opalDist);
45
51 void testNumEmitParticles(size_type nsteps, double dt) override;
52
58 void testEmitParticles(size_type nsteps, double dt) override;
59
60private:
61 using size_type = ippl::detail::size_type;
65 double distArea_m;
66 double sigmaTFall_m;
67 double sigmaTRise_m;
69 double fallTime_m;
70 double riseTime_m;
77
82 void setWithDomainDecomp(bool withDomainDecomp) override;
83
88 static size_t determineRandInit();
89
94 void setParameters(const std::shared_ptr<Distribution_t> &opalDist);
95
96public:
102 void generateUniformDisk(size_type nlocal, size_t nNew);
103
109
115 void generateParticles(size_t& numberOfParticles, Vector_t<double, 3> nr) override;
116
122 double FlatTopProfile(double t);
123
129 size_t computeNlocalUniformly(size_t nglobal);
130
139 double integrateTrapezoidal(double x1, double x2, double y1, double y2);
140
145 void initDomainDecomp(double BoxIncr) override;
146
153 double countEnteringParticlesPerRank(double t0, double tf);
154
159 void allocateParticles(size_t numberOfParticles);
160
166 void emitParticles(double t, double dt) override;
167};
168
169#endif // IPPL_FLAT_TOP_H
170
ParticleContainer< double, 3 > ParticleContainer_t
Definition Component.h:30
ippl::Vector< T, Dim > Vector_t
typename Kokkos::Random_XorShift64_Pool<> GeneratorPool
Distribution Distribution_t
Definition FlatTop.cpp:9
ippl::random::NormalDistribution< double, 3 > Dist_t
Definition FlatTop.cpp:11
FieldContainer< double, 3 > FieldContainer_t
const int nr
void setNr(Vector_t< double, 3 > nr)
Sets the number of grid points per direction.
Definition FlatTop.cpp:95
void testEmitParticles(size_type nsteps, double dt) override
Tests particle emission over a given number of steps.
Definition FlatTop.cpp:328
double distArea_m
Total area of the flattop distribution.
Definition FlatTop.h:65
double FlatTopProfile(double t)
Computes the flat-top profile value at a given time.
Definition FlatTop.cpp:116
bool emitting_m
Flag for particle emission status.
Definition FlatTop.h:71
Vector_t< double, 3 > nr_m
Number of grid points per direction.
Definition FlatTop.h:75
static size_t determineRandInit()
Determines the random seed initialization.
Definition FlatTop.cpp:24
void allocateParticles(size_t numberOfParticles)
Allocates memory for a given number of particles.
Definition FlatTop.cpp:254
double emissionTime_m
Total emission time.
Definition FlatTop.h:74
double sigmaTRise_m
Standard deviation for rise time profile.
Definition FlatTop.h:67
size_type totalN_m
Total number of particles.
Definition FlatTop.h:72
GeneratorPool rand_pool_m
Random number generator pool.
Definition FlatTop.h:62
double sigmaTFall_m
Standard deviation for fall time profile.
Definition FlatTop.h:66
bool withDomainDecomp_m
Flag for domain decomposition.
Definition FlatTop.h:73
void initDomainDecomp(double BoxIncr) override
Initializes the domain decomposition.
Definition FlatTop.cpp:157
Vector_t< double, 3 > hr_m
Grid spacing.
Definition FlatTop.h:76
double integrateTrapezoidal(double x1, double x2, double y1, double y2)
Integrates using the trapezoidal rule.
Definition FlatTop.cpp:153
double fallTime_m
Time duration for the fall phase.
Definition FlatTop.h:69
size_t computeNlocalUniformly(size_t nglobal)
Computes the local number of particles uniformly distributed among ranks.
Definition FlatTop.cpp:135
double normalizedFlankArea_m
Normalized area of the distribution flanks.
Definition FlatTop.h:64
void setWithDomainDecomp(bool withDomainDecomp) override
Sets whether to use domain decomposition.
Definition FlatTop.cpp:20
void emitParticles(double t, double dt) override
Emits new particles within a given time interval.
Definition FlatTop.cpp:264
double riseTime_m
Time duration for the rise phase.
Definition FlatTop.h:70
void generateParticles(size_t &numberOfParticles, Vector_t< double, 3 > nr) override
Generates particles with a given number and grid configuration.
Definition FlatTop.cpp:99
double countEnteringParticlesPerRank(double t0, double tf)
Counts the number of particles entering per rank in a given time interval.
Definition FlatTop.cpp:178
void generateUniformDisk(size_type nlocal, size_t nNew)
Generates particles (x,y) uniformly on a disk distribution.
Definition FlatTop.cpp:65
void testNumEmitParticles(size_type nsteps, double dt) override
Tests the number of emitted particles over a given number of steps.
Definition FlatTop.cpp:283
Vector_t< double, 3 > cutoffR_m
Cutoff radius.
Definition FlatTop.h:68
double flattopTime_m
Time duration of when the time profile is flat.
Definition FlatTop.h:63
FlatTop(std::shared_ptr< ParticleContainer_t > &pc, std::shared_ptr< FieldContainer_t > &fc, std::shared_ptr< Distribution_t > &opalDist)
Constructor for FlatTop.
Definition FlatTop.cpp:13
ippl::detail::size_type size_type
Definition FlatTop.h:61
void setParameters(const std::shared_ptr< Distribution_t > &opalDist)
Sets distribution parameters.
Definition FlatTop.cpp:36
SamplingBase(std::shared_ptr< ParticleContainer_t > &pc, std::shared_ptr< FieldContainer_t > &fc, std::shared_ptr< Distribution_t > &dist)