|
OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
|
Implements the sampling method for the flat-top distribution. x and y coordinates are uniformly distributed inside a circle and number of particles entering domain in [t, t+dt] follows flattop profile. More...
#include <FlatTop.h>
Public Member Functions | |
| FlatTop (std::shared_ptr< ParticleContainer_t > &pc, std::shared_ptr< FieldContainer_t > &fc, std::shared_ptr< Distribution_t > &opalDist) | |
| Constructor for FlatTop. | |
| void | testNumEmitParticles (size_type nsteps, double dt) override |
| Tests the number of emitted particles over a given number of steps. | |
| void | testEmitParticles (size_type nsteps, double dt) override |
| Tests particle emission over a given number of steps. | |
| void | generateUniformDisk (size_type nlocal, size_t nNew) |
| Generates particles (x,y) uniformly on a disk distribution. | |
| void | setNr (Vector_t< double, 3 > nr) |
| Sets the number of grid points per direction. | |
| void | generateParticles (size_t &numberOfParticles, Vector_t< double, 3 > nr) override |
| Generates particles with a given number and grid configuration. | |
| double | FlatTopProfile (double t) |
| Computes the flat-top profile value at a given time. | |
| size_t | computeNlocalUniformly (size_t nglobal) |
| Computes the local number of particles uniformly distributed among ranks. | |
| double | integrateTrapezoidal (double x1, double x2, double y1, double y2) |
| Integrates using the trapezoidal rule. | |
| void | initDomainDecomp (double BoxIncr) override |
| Initializes the domain decomposition. | |
| double | countEnteringParticlesPerRank (double t0, double tf) |
| Counts the number of particles entering per rank in a given time interval. | |
| void | allocateParticles (size_t numberOfParticles) |
| Allocates memory for a given number of particles. | |
| void | emitParticles (double t, double dt) override |
| Emits new particles within a given time interval. | |
| virtual void | testNumEmitParticles (size_t nsteps, double dt) |
| virtual void | testEmitParticles (size_t nsteps, double dt) |
Protected Attributes | |
| std::shared_ptr< ParticleContainer_t > | pc_m |
| std::shared_ptr< FieldContainer_t > | fc_m |
| std::shared_ptr< Distribution_t > | opalDist_m |
| std::string | samplingMethod_m |
Private Types | |
| using | size_type = ippl::detail::size_type |
Private Member Functions | |
| void | setWithDomainDecomp (bool withDomainDecomp) override |
| Sets whether to use domain decomposition. | |
| void | setParameters (const std::shared_ptr< Distribution_t > &opalDist) |
| Sets distribution parameters. | |
Static Private Member Functions | |
| static size_t | determineRandInit () |
| Determines the random seed initialization. | |
Private Attributes | |
| GeneratorPool | rand_pool_m |
| Random number generator pool. | |
| double | flattopTime_m |
| Time duration of when the time profile is flat. | |
| double | normalizedFlankArea_m |
| Normalized area of the distribution flanks. | |
| double | distArea_m |
| Total area of the flattop distribution. | |
| double | sigmaTFall_m |
| Standard deviation for fall time profile. | |
| double | sigmaTRise_m |
| Standard deviation for rise time profile. | |
| Vector_t< double, 3 > | cutoffR_m |
| Cutoff radius. | |
| double | fallTime_m |
| Time duration for the fall phase. | |
| double | riseTime_m |
| Time duration for the rise phase. | |
| bool | emitting_m |
| Flag for particle emission status. | |
| size_type | totalN_m |
| Total number of particles. | |
| bool | withDomainDecomp_m |
| Flag for domain decomposition. | |
| double | emissionTime_m |
| Total emission time. | |
| Vector_t< double, 3 > | nr_m |
| Number of grid points per direction. | |
| Vector_t< double, 3 > | hr_m |
| Grid spacing. | |
Implements the sampling method for the flat-top distribution. x and y coordinates are uniformly distributed inside a circle and number of particles entering domain in [t, t+dt] follows flattop profile.
The FlatTop distribution is f(t)/Z = exp[ -((t-riseTime_m)/sigma)^2/2 ] t < riseTime 1.0 riseTime < t < t<riseTime + flattopTime exp[ -((t-(fallTime_m + flattopTime_m))/sigmaTFall_m)^2/2 ] t>riseTime + flattopTime where Z is the normalizing factor.
|
private |
| FlatTop::FlatTop | ( | std::shared_ptr< ParticleContainer_t > & | pc, |
| std::shared_ptr< FieldContainer_t > & | fc, | ||
| std::shared_ptr< Distribution_t > & | opalDist ) |
Constructor for FlatTop.
| pc | Shared pointer to ParticleContainer. |
| fc | Shared pointer to FieldContainer. |
| opalDist | Shared pointer to Distribution. |
Definition at line 13 of file FlatTop.cpp.
References determineRandInit(), rand_pool_m, SamplingBase::SamplingBase(), and setParameters().
| void FlatTop::allocateParticles | ( | size_t | numberOfParticles | ) |
Allocates memory for a given number of particles.
| numberOfParticles | Number of particles to allocate. |
Definition at line 254 of file FlatTop.cpp.
References computeNlocalUniformly(), SamplingBase::pc_m, and totalN_m.
Referenced by generateParticles().
| size_t FlatTop::computeNlocalUniformly | ( | size_t | nglobal | ) |
Computes the local number of particles uniformly distributed among ranks.
| nglobal | Total global number of particles. |
Definition at line 135 of file FlatTop.cpp.
Referenced by allocateParticles(), and countEnteringParticlesPerRank().
| double FlatTop::countEnteringParticlesPerRank | ( | double | t0, |
| double | tf ) |
Counts the number of particles entering per rank in a given time interval.
| t0 | Start time. |
| tf | End time. |
Definition at line 178 of file FlatTop.cpp.
References Physics::c, computeNlocalUniformly(), Dim, distArea_m, FlatTopProfile(), integrateTrapezoidal(), SamplingBase::opalDist_m, SamplingBase::pc_m, totalN_m, and withDomainDecomp_m.
Referenced by emitParticles(), and testNumEmitParticles().
|
staticprivate |
Determines the random seed initialization.
Definition at line 24 of file FlatTop.cpp.
References gmsg, and Options::seed.
Referenced by FlatTop().
|
overridevirtual |
Emits new particles within a given time interval.
| t | Start time. |
| dt | Time step. |
Reimplemented from SamplingBase.
Definition at line 264 of file FlatTop.cpp.
References countEnteringParticlesPerRank(), generateUniformDisk(), gmsg, and SamplingBase::pc_m.
Referenced by testEmitParticles().
| double FlatTop::FlatTopProfile | ( | double | t | ) |
Computes the flat-top profile value at a given time.
| t | Time value. |
Definition at line 116 of file FlatTop.cpp.
References fallTime_m, flattopTime_m, riseTime_m, sigmaTFall_m, and sigmaTRise_m.
Referenced by countEnteringParticlesPerRank().
|
overridevirtual |
Generates particles with a given number and grid configuration.
| numberOfParticles | Number of particles to generate. |
| nr | Number of grid points in each dimension. |
Reimplemented from SamplingBase.
Definition at line 99 of file FlatTop.cpp.
References allocateParticles(), emitting_m, nr, SamplingBase::pc_m, and setNr().
| void FlatTop::generateUniformDisk | ( | size_type | nlocal, |
| size_t | nNew ) |
Generates particles (x,y) uniformly on a disk distribution.
| nlocal | Number of local particles. |
| nNew | Number of new particles to generate. |
Definition at line 65 of file FlatTop.cpp.
References SamplingBase::opalDist_m, SamplingBase::pc_m, Physics::pi, pi, and rand_pool_m.
Referenced by emitParticles(), and testNumEmitParticles().
|
overridevirtual |
Initializes the domain decomposition.
| BoxIncr | Box increment factor. |
Reimplemented from SamplingBase.
Definition at line 157 of file FlatTop.cpp.
References Physics::c, emissionTime_m, SamplingBase::fc_m, hr_m, nr_m, SamplingBase::opalDist_m, and SamplingBase::pc_m.
| double FlatTop::integrateTrapezoidal | ( | double | x1, |
| double | x2, | ||
| double | y1, | ||
| double | y2 ) |
Integrates using the trapezoidal rule.
| x1 | begining of interval [x1,x2]. |
| x2 | end of interval [x1,x2]. |
| y1 | value of f(x1). |
| y2 | value of f(x2). |
Definition at line 153 of file FlatTop.cpp.
Referenced by countEnteringParticlesPerRank().
| void FlatTop::setNr | ( | Vector_t< double, 3 > | nr | ) |
Sets the number of grid points per direction.
| nr | Vector specifying the number of grid points. |
Definition at line 95 of file FlatTop.cpp.
Referenced by generateParticles().
|
private |
Sets distribution parameters.
| opalDist | Shared pointer to the distribution object. |
Definition at line 36 of file FlatTop.cpp.
References cutoffR_m, distArea_m, emissionTime_m, emitting_m, fallTime_m, SamplingBase::fc_m, flattopTime_m, normalizedFlankArea_m, SamplingBase::opalDist_m, riseTime_m, sigmaTFall_m, sigmaTRise_m, and Physics::two_pi.
Referenced by FlatTop().
|
overrideprivatevirtual |
Sets whether to use domain decomposition.
| withDomainDecomp | Boolean flag for domain decomposition. |
Reimplemented from SamplingBase.
Definition at line 20 of file FlatTop.cpp.
References withDomainDecomp_m.
|
override |
Tests particle emission over a given number of steps.
| nsteps | Number of time steps to simulate. |
| dt | Time step size. |
Definition at line 328 of file FlatTop.cpp.
References emitParticles().
|
inlinevirtualinherited |
Definition at line 37 of file SamplingBase.hpp.
|
override |
Tests the number of emitted particles over a given number of steps.
| nsteps | Number of time steps to simulate. |
| dt | Time step size. |
Definition at line 283 of file FlatTop.cpp.
References Physics::c, countEnteringParticlesPerRank(), emissionTime_m, generateUniformDisk(), and SamplingBase::pc_m.
|
inlinevirtualinherited |
Definition at line 34 of file SamplingBase.hpp.
|
private |
|
private |
Total area of the flattop distribution.
Definition at line 65 of file FlatTop.h.
Referenced by countEnteringParticlesPerRank(), and setParameters().
|
private |
Total emission time.
Definition at line 74 of file FlatTop.h.
Referenced by initDomainDecomp(), setParameters(), and testNumEmitParticles().
|
private |
Flag for particle emission status.
Definition at line 71 of file FlatTop.h.
Referenced by generateParticles(), and setParameters().
|
private |
Time duration for the fall phase.
Definition at line 69 of file FlatTop.h.
Referenced by FlatTopProfile(), and setParameters().
|
protectedinherited |
Definition at line 14 of file SamplingBase.hpp.
Referenced by FlatTop::initDomainDecomp(), SamplingBase(), and FlatTop::setParameters().
|
private |
Time duration of when the time profile is flat.
Definition at line 63 of file FlatTop.h.
Referenced by FlatTopProfile(), and setParameters().
|
private |
|
private |
Normalized area of the distribution flanks.
Definition at line 64 of file FlatTop.h.
Referenced by setParameters().
|
private |
Number of grid points per direction.
Definition at line 75 of file FlatTop.h.
Referenced by initDomainDecomp(), and setNr().
|
protectedinherited |
Definition at line 15 of file SamplingBase.hpp.
Referenced by FlatTop::countEnteringParticlesPerRank(), FlatTop::generateUniformDisk(), FlatTop::initDomainDecomp(), MultiVariateGaussian::MultiVariateGaussian(), SamplingBase(), and FlatTop::setParameters().
|
protectedinherited |
Definition at line 13 of file SamplingBase.hpp.
Referenced by FlatTop::allocateParticles(), FlatTop::countEnteringParticlesPerRank(), FlatTop::emitParticles(), FlatTop::generateParticles(), Gaussian::generateParticles(), MultiVariateGaussian::generateParticles(), FlatTop::generateUniformDisk(), FlatTop::initDomainDecomp(), SamplingBase(), SamplingBase(), and FlatTop::testNumEmitParticles().
|
private |
Random number generator pool.
Definition at line 62 of file FlatTop.h.
Referenced by FlatTop(), and generateUniformDisk().
|
private |
Time duration for the rise phase.
Definition at line 70 of file FlatTop.h.
Referenced by FlatTopProfile(), and setParameters().
|
protectedinherited |
Definition at line 16 of file SamplingBase.hpp.
|
private |
Standard deviation for fall time profile.
Definition at line 66 of file FlatTop.h.
Referenced by FlatTopProfile(), and setParameters().
|
private |
Standard deviation for rise time profile.
Definition at line 67 of file FlatTop.h.
Referenced by FlatTopProfile(), and setParameters().
|
private |
Total number of particles.
Definition at line 72 of file FlatTop.h.
Referenced by allocateParticles(), and countEnteringParticlesPerRank().
|
private |
Flag for domain decomposition.
Definition at line 73 of file FlatTop.h.
Referenced by countEnteringParticlesPerRank(), and setWithDomainDecomp().