OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
FlatTop Class Reference

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>

Inheritance diagram for FlatTop:
Inheritance graph
Collaboration diagram for FlatTop:
Collaboration graph

Public Member Functions

 FlatTop (std::shared_ptr< ParticleContainer_t > pc, std::shared_ptr< FieldContainer_t > fc, Distribution_t *opalDist)
 Constructor for FlatTop.
 
 FlatTop (std::shared_ptr< ParticleContainer_t > pc, std::shared_ptr< FieldContainer_t > fc, bool emitting, double sigmaTFall, double sigmaTRise, Vector_t< double, 3 > cutoff, double tPulseLengthFWHM, Vector_t< double, 3 > sigmaR)
 Constructor for FlatTop.
 
 FlatTop (std::shared_ptr< ParticleContainer_t > pc, bool emitting, double sigmaTFall, double sigmaTRise, Vector_t< double, 3 > cutoff, double tPulseLengthFWHM, Vector_t< double, 3 > sigmaR)
 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 setWithDomainDecomp (bool withDomainDecomp) override
 Sets whether to use domain decomposition.
 
void setDistArea (double distArea)
 Sets the total area of the distribution.
 
double getDistArea ()
 Returns the total area of the distribution.
 
double getEmissionTime ()
 Returns the total emission time.
 
bool isEmissionDone (double t) const override
 Whether this sampler has finished all emission (no more particles will be created).
 
void generateUniformDisk (size_type nlocal, size_t nNew, double dt)
 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.
 
size_type 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.
 
void setInternalVariables (bool emitting, double sigmaTFall, double sigmaTRise, Vector_t< double, 3 > cutoff, double tPulseLengthFWHM, Vector_t< double, 3 > sigmaR)
 
void setEmissionOffsets (ippl::Vector< double, 3 > R0, ippl::Vector< double, 3 > P0, double t0, const std::string &emissionModel="NONE")
 
void setInitialPolarization (const ippl::Vector< double, 3 > &pol)
 
void fillPolarization (size_t startIdx, size_t count)
 
Vector_t< double, 3 > getEmissionR0 () const
 
virtual double getGlobalTimeShift () const
 Optional tracker time shift needed before the first emitted particles are born.
 
virtual bool hasInitialReferenceMomentum () const
 Optional initial reference momentum in the source-local frame.
 
virtual Vector_t< double, 3 > getInitialReferenceMomentum () const
 
virtual double getEmissionTimeStep () const
 Optional OPAL-like tracker time step while this sampler is still emitting.
 
virtual void testNumEmitParticles (size_t, double)
 
virtual void testEmitParticles (size_t, double)
 
size_t computeLocalEmitCount (size_t totalToSample) const
 Computes the number of particles this rank should emit so that the global total equals totalToSample and no rank exceeds its capacity (space left).
 

Protected Attributes

std::shared_ptr< ParticleContainer_tpc_m
 
std::shared_ptr< FieldContainer_tfc_m
 
Distribution_topalDist_m = nullptr
 
std::string samplingMethod_m
 
Vector_t< double, 3 > R0_m = 0.0
 Emission source offset: position R0, momentum P0, start time t0 (applied in sample step).
 
Vector_t< double, 3 > P0_m = 0.0
 
double t0_m = 0.0
 
std::string emissionModel_m = "NONE"
 
Vector_t< double, 3 > initialPol_m = 0.0
 
bool hasEmittedOnce_m = false
 For one-shot emitters (e.g. Gaussian at delayed t0): guard to avoid double sampling.
 

Private Types

using size_type = ippl::detail::size_type
 

Private Member Functions

void setParameters (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.
 
Vector_t< double, 3 > sigmaR_m
 

Detailed Description

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.

Definition at line 33 of file FlatTop.h.

Member Typedef Documentation

◆ size_type

using FlatTop::size_type = ippl::detail::size_type
private

Definition at line 114 of file FlatTop.h.

Constructor & Destructor Documentation

◆ FlatTop() [1/3]

FlatTop::FlatTop ( std::shared_ptr< ParticleContainer_t pc,
std::shared_ptr< FieldContainer_t fc,
Distribution_t opalDist 
)

Constructor for FlatTop.

Parameters
pcShared pointer to ParticleContainer.
fcShared pointer to FieldContainer.
opalDistBorrowed Distribution.

Definition at line 10 of file FlatTop.cpp.

References setParameters().

Here is the call graph for this function:

◆ FlatTop() [2/3]

FlatTop::FlatTop ( std::shared_ptr< ParticleContainer_t pc,
std::shared_ptr< FieldContainer_t fc,
bool  emitting,
double  sigmaTFall,
double  sigmaTRise,
Vector_t< double, 3 >  cutoff,
double  tPulseLengthFWHM,
Vector_t< double, 3 >  sigmaR 
)

Constructor for FlatTop.

Parameters
pcShared pointer to ParticleContainer.
fcShared pointer to FieldContainer.
emittingFlag indicating whether the distribution is emitting.
sigmaTFallStandard deviation of fall in flat-top profile.
sigmaTRiseStandard deviation of rise in flat-top profile.
cutoffCutoff multiplier for position distribution R.
tPulseLengthFWHMTime length of the pulse in FWHM.
sigmaRStandard deviation for position distribution.

Definition at line 17 of file FlatTop.cpp.

References setInternalVariables().

Here is the call graph for this function:

◆ FlatTop() [3/3]

FlatTop::FlatTop ( std::shared_ptr< ParticleContainer_t pc,
bool  emitting,
double  sigmaTFall,
double  sigmaTRise,
Vector_t< double, 3 >  cutoff,
double  tPulseLengthFWHM,
Vector_t< double, 3 >  sigmaR 
)

Constructor for FlatTop.

Parameters
pcShared pointer to ParticleContainer.
emittingFlag indicating whether the distribution is emitting.
sigmaTFallStandard deviation of fall in flat-top profile.
sigmaTRiseStandard deviation of rise in flat-top profile.
cutoffCutoff multiplier for position distribution R.
tPulseLengthFWHMTime span of the flat-top profile.
sigmaRStandard deviation for position distribution R.

Definition at line 25 of file FlatTop.cpp.

References setInternalVariables().

Here is the call graph for this function:

Member Function Documentation

◆ allocateParticles()

void FlatTop::allocateParticles ( size_t  numberOfParticles)

Allocates memory for a given number of particles.

Parameters
numberOfParticlesNumber of particles to allocate.

Definition at line 350 of file FlatTop.cpp.

References totalN_m.

Referenced by generateParticles(), and TEST_F().

◆ computeLocalEmitCount()

size_t SamplingBase::computeLocalEmitCount ( size_t  totalToSample) const
inherited

Computes the number of particles this rank should emit so that the global total equals totalToSample and no rank exceeds its capacity (space left).

Parameters
totalToSampleGlobal number of particles to emit this timestep.
Returns
Local number of particles this rank should create/emit.

Definition at line 6 of file SamplingBase.cpp.

References SamplingBase::pc_m.

Referenced by countEnteringParticlesPerRank(), FromFile::generateParticles(), Gaussian::generateParticles(), and MultiVariateGaussian::generateParticles().

◆ computeNlocalUniformly()

size_t FlatTop::computeNlocalUniformly ( size_t  nglobal)

Computes the local number of particles uniformly distributed among ranks.

Parameters
nglobalTotal global number of particles.
Returns
Number of local particles.

Definition at line 222 of file FlatTop.cpp.

Referenced by countEnteringParticlesPerRank().

◆ countEnteringParticlesPerRank()

FlatTop::size_type FlatTop::countEnteringParticlesPerRank ( double  t0,
double  tf 
)

Counts the number of particles entering per rank in a given time interval.

Parameters
t0Start time.
tfEnd time.
Returns
Number of entering particles per rank.

Definition at line 263 of file FlatTop.cpp.

References SamplingBase::computeLocalEmitCount(), computeNlocalUniformly(), distArea_m, FlatTopProfile(), integrateTrapezoidal(), SamplingBase::pc_m, and totalN_m.

Referenced by emitParticles(), TEST_F(), and testNumEmitParticles().

Here is the call graph for this function:

◆ determineRandInit()

size_t FlatTop::determineRandInit ( )
staticprivate

Determines the random seed initialization.

Returns
The seed value.

Definition at line 35 of file FlatTop.cpp.

References gmsg, and Options::seed.

◆ emitParticles()

void FlatTop::emitParticles ( double  t,
double  dt 
)
overridevirtual

Emits new particles within a given time interval.

Parameters
tStart time.
dtTime step.

Reimplemented from SamplingBase.

Definition at line 360 of file FlatTop.cpp.

References countEnteringParticlesPerRank(), SamplingBase::fillPolarization(), generateUniformDisk(), SamplingBase::pc_m, and SamplingBase::t0_m.

Referenced by testEmitParticles().

Here is the call graph for this function:

◆ fillPolarization()

void SamplingBase::fillPolarization ( size_t  startIdx,
size_t  count 
)
inlineinherited

Fill the Pol particle attribute with initialPol_m on the half-open range [startIdx, startIdx + count). No-op when the container does not store spin. Must be called after the particles have been created and before they are pushed downstream.

Definition at line 63 of file SamplingBase.hpp.

References SamplingBase::initialPol_m, and SamplingBase::pc_m.

Referenced by emitParticles(), FromFile::generateParticles(), Gaussian::generateParticles(), and MultiVariateGaussian::generateParticles().

◆ FlatTopProfile()

double FlatTop::FlatTopProfile ( double  t)

Computes the flat-top profile value at a given time.

Parameters
tTime value.
Returns
Profile value at time t.

Definition at line 204 of file FlatTop.cpp.

References fallTime_m, flattopTime_m, riseTime_m, sigmaTFall_m, and sigmaTRise_m.

Referenced by countEnteringParticlesPerRank(), and TEST_F().

◆ generateParticles()

void FlatTop::generateParticles ( size_t &  numberOfParticles,
Vector_t< double, 3 >  nr 
)
overridevirtual

Generates particles with a given number and grid configuration.

Parameters
numberOfParticlesNumber of particles to generate.
nrNumber of grid points in each dimension.

Reimplemented from SamplingBase.

Definition at line 183 of file FlatTop.cpp.

References allocateParticles(), nr, and setNr().

Here is the call graph for this function:

◆ generateUniformDisk()

void FlatTop::generateUniformDisk ( size_type  nlocal,
size_t  nNew,
double  dt 
)

Generates particles (x,y) uniformly on a disk distribution.

Each new particle is assigned a fractional per-particle dt sampled uniformly in (0, dt), representing the fraction of the next integration step the particle will experience (as if born at a random time within the current emission interval).

Parameters
nlocalNumber of local particles.
nNewNumber of new particles to generate.
dtGlobal timestep; used to sample fractional per-particle dt.

Definition at line 88 of file FlatTop.cpp.

References Physics::c, SamplingBase::emissionModel_m, euclidean_norm(), SamplingBase::P0_m, SamplingBase::pc_m, Physics::pi, SamplingBase::R0_m, rand_pool_m, and sigmaR_m.

Referenced by emitParticles(), TEST_F(), and testNumEmitParticles().

Here is the call graph for this function:

◆ getDistArea()

double FlatTop::getDistArea ( )
inline

Returns the total area of the distribution.

Definition at line 103 of file FlatTop.h.

References distArea_m.

Referenced by TEST_F().

◆ getEmissionR0()

Vector_t< double, 3 > SamplingBase::getEmissionR0 ( ) const
inlineinherited

Definition at line 78 of file SamplingBase.hpp.

References SamplingBase::R0_m.

◆ getEmissionTime()

double FlatTop::getEmissionTime ( )
inline

Returns the total emission time.

Definition at line 108 of file FlatTop.h.

References emissionTime_m.

◆ getEmissionTimeStep()

virtual double SamplingBase::getEmissionTimeStep ( ) const
inlinevirtualinherited

Optional OPAL-like tracker time step while this sampler is still emitting.

Reimplemented in EmittedFromFile, and OpalFlatTop.

Definition at line 102 of file SamplingBase.hpp.

◆ getGlobalTimeShift()

virtual double SamplingBase::getGlobalTimeShift ( ) const
inlinevirtualinherited

Optional tracker time shift needed before the first emitted particles are born.

Reimplemented in EmittedFromFile, and OpalFlatTop.

Definition at line 92 of file SamplingBase.hpp.

◆ getInitialReferenceMomentum()

virtual Vector_t< double, 3 > SamplingBase::getInitialReferenceMomentum ( ) const
inlinevirtualinherited

Reimplemented in EmittedFromFile, and OpalFlatTop.

Definition at line 97 of file SamplingBase.hpp.

◆ hasInitialReferenceMomentum()

virtual bool SamplingBase::hasInitialReferenceMomentum ( ) const
inlinevirtualinherited

Optional initial reference momentum in the source-local frame.

Reimplemented in EmittedFromFile, and OpalFlatTop.

Definition at line 95 of file SamplingBase.hpp.

◆ initDomainDecomp()

void FlatTop::initDomainDecomp ( double  BoxIncr)
overridevirtual

Initializes the domain decomposition.

Parameters
BoxIncrBox increment factor.

Reimplemented from SamplingBase.

Definition at line 242 of file FlatTop.cpp.

References Physics::c, emissionTime_m, SamplingBase::fc_m, hr_m, nr_m, SamplingBase::pc_m, and sigmaR_m.

◆ integrateTrapezoidal()

double FlatTop::integrateTrapezoidal ( double  x1,
double  x2,
double  y1,
double  y2 
)

Integrates using the trapezoidal rule.

Parameters
x1begining of interval [x1,x2].
x2end of interval [x1,x2].
y1value of f(x1).
y2value of f(x2).
Returns
Integrated result.

Definition at line 238 of file FlatTop.cpp.

Referenced by countEnteringParticlesPerRank().

◆ isEmissionDone()

bool FlatTop::isEmissionDone ( double  t) const
inlineoverridevirtual

Whether this sampler has finished all emission (no more particles will be created).

Parameters
tCurrent simulation time (s).

Reimplemented from SamplingBase.

Definition at line 111 of file FlatTop.h.

References emissionTime_m, and SamplingBase::t0_m.

◆ setDistArea()

void FlatTop::setDistArea ( double  distArea)
inline

Sets the total area of the distribution.

Definition at line 97 of file FlatTop.h.

References distArea_m.

◆ setEmissionOffsets()

void SamplingBase::setEmissionOffsets ( ippl::Vector< double, 3 >  R0,
ippl::Vector< double, 3 >  P0,
double  t0,
const std::string &  emissionModel = "NONE" 
)
inlineinherited

◆ setInitialPolarization()

void SamplingBase::setInitialPolarization ( const ippl::Vector< double, 3 > &  pol)
inlineinherited

Definition at line 57 of file SamplingBase.hpp.

References SamplingBase::initialPol_m.

◆ setInternalVariables()

void FlatTop::setInternalVariables ( bool  emitting,
double  sigmaTFall,
double  sigmaTRise,
Vector_t< double, 3 >  cutoff,
double  tPulseLengthFWHM,
Vector_t< double, 3 >  sigmaR 
)

◆ setNr()

void FlatTop::setNr ( Vector_t< double, 3 >  nr)

Sets the number of grid points per direction.

Parameters
nrVector specifying the number of grid points.

Definition at line 181 of file FlatTop.cpp.

References nr, and nr_m.

Referenced by generateParticles().

◆ setParameters()

void FlatTop::setParameters ( Distribution_t opalDist)
private

Sets distribution parameters.

Parameters
opalDistShared pointer to the distribution object.

Definition at line 47 of file FlatTop.cpp.

References emissionTime_m, Distribution::emitting_m, SamplingBase::fc_m, Distribution::getCutoffR(), Distribution::getSigmaR(), Distribution::getSigmaTFall(), Distribution::getSigmaTRise(), Distribution::getTPulseLengthFWHM(), SamplingBase::opalDist_m, setInternalVariables(), and Distribution::setTEmission().

Referenced by FlatTop().

Here is the call graph for this function:

◆ setWithDomainDecomp()

void FlatTop::setWithDomainDecomp ( bool  withDomainDecomp)
overridevirtual

Sets whether to use domain decomposition.

Parameters
withDomainDecompBoolean flag for domain decomposition.

Reimplemented from SamplingBase.

Definition at line 33 of file FlatTop.cpp.

References withDomainDecomp_m.

Referenced by TEST_F().

◆ testEmitParticles() [1/2]

virtual void SamplingBase::testEmitParticles ( size_t  ,
double   
)
inlinevirtualinherited

Definition at line 108 of file SamplingBase.hpp.

◆ testEmitParticles() [2/2]

void FlatTop::testEmitParticles ( size_type  nsteps,
double  dt 
)
override

Tests particle emission over a given number of steps.

Parameters
nstepsNumber of time steps to simulate.
dtTime step size.

Definition at line 440 of file FlatTop.cpp.

References emitParticles().

Here is the call graph for this function:

◆ testNumEmitParticles() [1/2]

virtual void SamplingBase::testNumEmitParticles ( size_t  ,
double   
)
inlinevirtualinherited

Definition at line 105 of file SamplingBase.hpp.

◆ testNumEmitParticles() [2/2]

void FlatTop::testNumEmitParticles ( size_type  nsteps,
double  dt 
)
override

Tests the number of emitted particles over a given number of steps.

Parameters
nstepsNumber of time steps to simulate.
dtTime step size.

Definition at line 394 of file FlatTop.cpp.

References Physics::c, countEnteringParticlesPerRank(), emissionTime_m, generateUniformDisk(), and SamplingBase::pc_m.

Here is the call graph for this function:

Member Data Documentation

◆ cutoffR_m

Vector_t<double, 3> FlatTop::cutoffR_m
private

Cutoff radius.

Definition at line 121 of file FlatTop.h.

Referenced by setInternalVariables().

◆ distArea_m

double FlatTop::distArea_m
private

Total area of the flattop distribution.

Definition at line 118 of file FlatTop.h.

Referenced by countEnteringParticlesPerRank(), getDistArea(), setDistArea(), and setInternalVariables().

◆ emissionModel_m

◆ emissionTime_m

double FlatTop::emissionTime_m
private

Total emission time.

Definition at line 127 of file FlatTop.h.

Referenced by getEmissionTime(), initDomainDecomp(), isEmissionDone(), setInternalVariables(), setParameters(), and testNumEmitParticles().

◆ emitting_m

bool FlatTop::emitting_m
private

Flag for particle emission status.

Definition at line 124 of file FlatTop.h.

Referenced by setInternalVariables().

◆ fallTime_m

double FlatTop::fallTime_m
private

Time duration for the fall phase.

Definition at line 122 of file FlatTop.h.

Referenced by FlatTopProfile(), and setInternalVariables().

◆ fc_m

std::shared_ptr<FieldContainer_t> SamplingBase::fc_m
protectedinherited

◆ flattopTime_m

double FlatTop::flattopTime_m
private

Time duration of when the time profile is flat.

Definition at line 116 of file FlatTop.h.

Referenced by FlatTopProfile(), and setInternalVariables().

◆ hasEmittedOnce_m

bool SamplingBase::hasEmittedOnce_m = false
protectedinherited

◆ hr_m

Vector_t<double, 3> FlatTop::hr_m
private

Grid spacing.

Definition at line 129 of file FlatTop.h.

Referenced by initDomainDecomp().

◆ initialPol_m

Vector_t<double, 3> SamplingBase::initialPol_m = 0.0
protectedinherited

Initial polarization vector P applied to every particle this sampler creates. Rest-frame Pauli expectations along lab-frame axes; |Pol| in [0, 1]. Ignored when the container does not have a spin attribute (hasSpin() == false).

Definition at line 30 of file SamplingBase.hpp.

Referenced by SamplingBase::fillPolarization(), and SamplingBase::setInitialPolarization().

◆ normalizedFlankArea_m

double FlatTop::normalizedFlankArea_m
private

Normalized area of the distribution flanks.

Definition at line 117 of file FlatTop.h.

Referenced by setInternalVariables().

◆ nr_m

Vector_t<double, 3> FlatTop::nr_m
private

Number of grid points per direction.

Definition at line 128 of file FlatTop.h.

Referenced by initDomainDecomp(), and setNr().

◆ opalDist_m

◆ P0_m

◆ pc_m

◆ R0_m

Vector_t<double, 3> SamplingBase::R0_m = 0.0
protectedinherited

◆ rand_pool_m

GeneratorPool FlatTop::rand_pool_m
private

Random number generator pool.

Definition at line 115 of file FlatTop.h.

Referenced by generateUniformDisk().

◆ riseTime_m

double FlatTop::riseTime_m
private

Time duration for the rise phase.

Definition at line 123 of file FlatTop.h.

Referenced by FlatTopProfile(), and setInternalVariables().

◆ samplingMethod_m

std::string SamplingBase::samplingMethod_m
protectedinherited

Definition at line 20 of file SamplingBase.hpp.

◆ sigmaR_m

Vector_t<double, 3> FlatTop::sigmaR_m
private

Definition at line 130 of file FlatTop.h.

Referenced by generateUniformDisk(), initDomainDecomp(), and setInternalVariables().

◆ sigmaTFall_m

double FlatTop::sigmaTFall_m
private

Standard deviation for fall time profile.

Definition at line 119 of file FlatTop.h.

Referenced by FlatTopProfile(), and setInternalVariables().

◆ sigmaTRise_m

double FlatTop::sigmaTRise_m
private

Standard deviation for rise time profile.

Definition at line 120 of file FlatTop.h.

Referenced by FlatTopProfile(), and setInternalVariables().

◆ t0_m

◆ totalN_m

size_type FlatTop::totalN_m
private

Total number of particles.

Definition at line 125 of file FlatTop.h.

Referenced by allocateParticles(), and countEnteringParticlesPerRank().

◆ withDomainDecomp_m

bool FlatTop::withDomainDecomp_m
private

Flag for domain decomposition.

Definition at line 126 of file FlatTop.h.

Referenced by setWithDomainDecomp().


The documentation for this class was generated from the following files: