milupHPC documentation
Public Member Functions | Public Attributes | Private Member Functions | List of all members
Leapfrog Class Reference

#include "leapfrog.h"

+ Inheritance diagram for Leapfrog:
Inheritance graph
+ Collaboration diagram for Leapfrog:
Collaboration graph

Public Member Functions

 Leapfrog (SimulationParameters simulationParameters)
 Constructor. More...
 
 ~Leapfrog ()
 Destructor. More...
 
void integrate (int step)
 Implementation of the abstract integration method. More...
 
real removeParticles ()
 Remove particles in dependence of some criterion. More...
 
void fixedLoadBalancing ()
 Load balancing via equidistant ranges. More...
 
void dynamicLoadBalancing (int bins=5000)
 Pre-calculations for updateRangeApproximately. More...
 
void updateRangeApproximately (int aimedParticlesPerProcess, int bins=5000)
 Update the ranges (approximately and dynamically). More...
 
void updateRange (int aimedParticlesPerProcess)
 Update the range in dependence on number of (MPI) processes and aimed particles per process. More...
 
void prepareSimulation ()
 Prepare the simulation, including. More...
 
void numParticlesFromFile (const std::string &filename)
 Determine amount of particles (numParticles and numParticlesLocal) from initial file/particle distribution file. More...
 
void distributionFromFile (const std::string &filename)
 Read initial/particle distribution file (in parallel) More...
 
real tree ()
 Wrapper function for building the tree (and domain tree). More...
 
real pseudoParticles ()
 Wrapper function for calculating pseudo-particles. More...
 
real gravity ()
 Wrapper function for Gravity-related stuff. More...
 
real sph ()
 
real rhs (int step, bool selfGravity=true, bool assignParticlesToProcess=true)
 
void afterIntegrationStep ()
 
real particles2file (int step)
 
real particles2file (const std::string &filename, int *particleIndices, int length)
 
real configParameters2file (HighFive::File &h5file)
 
void getMemoryInfo ()
 

Public Attributes

int subStep
 current sub-step (there are possibly more sub-steps within a step!) More...
 
real h_searchRadius
 search radius for SPH (MPI-process overarching) neighbor search More...
 
real totalEnergy
 total energy More...
 
real totalAngularMomentum
 total angular momentum More...
 
H5Profiler & profiler = H5Profiler::getInstance("log/performance.h5")
 H5 profiler instance. More...
 
Curve::Type curveType
 Space-filling curve type to be used (Lebesgue or Hilbert) More...
 
SPH::KernelHandler kernelHandler
 Instance to handle the SPH Kernel instance on device and host. More...
 
SimulationTimeHandler * simulationTimeHandler
 Instance to handle the SimulationTime instances on device and host. More...
 
integer numParticles
 number of particles (to be allocated) More...
 
integer sumParticles
 (real) number of particles on all processes More...
 
integer numParticlesLocal
 
integer numNodes
 number of nodes (to be allocated) More...
 
IntegratedParticleHandler * integratedParticles
 
integer * d_mutex
 
ParticleHandler * particleHandler
 Instance to handle the Particles instance on device and host. More...
 
SubDomainKeyTreeHandler * subDomainKeyTreeHandler
 Instance to handle the SubDomainKeyTree instance on device and host. More...
 
TreeHandler * treeHandler
 Instance to handle the Tree instance on device and host. More...
 
DomainListHandler * domainListHandler
 Instance to handle the DomainList instance on device and host. More...
 
DomainListHandler * lowestDomainListHandler
 Instance to handle the (lowest) DomainList instance on device and host. More...
 
MaterialHandler * materialHandler
 Instance to handle Materials instances on device and host. More...
 
HelperHandler * buffer
 buffer instance More...
 
SimulationParameters simulationParameters
 buffer (need for revising) More...
 

Private Member Functions

real reset ()
 Reset arrays, values, ... More...
 
real boundingBox ()
 Calculate bounding boxes/simulation domain. More...
 
real parallel_tree ()
 Parallel version regarding tree-stuff. More...
 
real parallel_pseudoParticles ()
 Parallel version regarding computation of pseudo-particles. More...
 
real parallel_gravity ()
 Parallel version regarding computation of gravitational stuff. More...
 
real parallel_sph ()
 Parallel version regarding computation of SPH-stuff. More...
 
template<typename T >
integer sendParticlesEntry (integer *sendLengths, integer *receiveLengths, T *entry, T *entryBuffer, T *copyBuffer)
 Send particles/Exchange particles among MPI processes. More...
 
template<typename T >
integer sendParticles (T *sendBuffer, T *receiveBuffer, integer *sendLengths, integer *receiveLengths)
 Send particles/Exchange particles among MPI processes. More...
 
template<typename U , typename T >
real arrangeParticleEntries (U *sortArray, U *sortedArray, T *entry, T *temp)
 Function to sort an array entry in dependence of another array sortArray More...
 
real assignParticles ()
 Assign particles to correct process in dependence of particle key and ranges. More...
 
real angularMomentum ()
 Calculate the angular momentum for all particles. More...
 
real energy ()
 Calculate the total amount of energy. More...
 

Detailed Description

Definition at line 16 of file leapfrog.h.

Constructor & Destructor Documentation

◆ Leapfrog()

Leapfrog::Leapfrog ( SimulationParameters  simulationParameters)

Constructor.

Parameters
simulationParametersSimulation parameters/settings.

Definition at line 9 of file leapfrog.cpp.

+ Here is the call graph for this function:

◆ ~Leapfrog()

Leapfrog::~Leapfrog ( )

Destructor.

Definition at line 15 of file leapfrog.cpp.

+ Here is the call graph for this function:

Member Function Documentation

◆ afterIntegrationStep()

void Miluphpc::afterIntegrationStep ( )
inherited

Function to be called after successful integration (step)

Definition at line 471 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ angularMomentum()

real Miluphpc::angularMomentum ( )
privateinherited

Calculate the angular momentum for all particles.

Returns
accumulated time of functions within

Definition at line 484 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ arrangeParticleEntries()

template<typename U , typename T >
real Miluphpc::arrangeParticleEntries ( U *  sortArray,
U *  sortedArray,
T *  entry,
T *  temp 
)
privateinherited

Function to sort an array entry in dependence of another array sortArray

Template Parameters
Utype of the arrays determining the sorting
Ttype of the arrays to be sorted
Parameters
sortArrayarray to be sorted (rather sorting behaviour should be determined)
sortedArraysorted array (result of sorting)
entry(relevant) array to be sorted
tempbuffer needed for sorting process
Returns
accumulated time of functions within

Definition at line 944 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ assignParticles()

real Miluphpc::assignParticles ( )
privateinherited

Assign particles to correct process in dependence of particle key and ranges.

First an extra array is used to save the information of process assignment for each particle, this extra array is used as key for sorting all of the attributes of the Particle class and finally the sub-array determined by the process assignment are sent to the corresponding process via MPI.

Returns
accumulated time of functions within

Definition at line 700 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ boundingBox()

real Miluphpc::boundingBox ( )
privateinherited

Calculate bounding boxes/simulation domain.

To derive the bounding boxes of the simulation domain the maximum and minimum value for each coordinate axis is determined locally on each process and afterwards reduced to all processes to get the global maximum and minimum via MPI_Allreduce().

Returns
accumulated time of functions within

Definition at line 643 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ configParameters2file()

real Miluphpc::configParameters2file ( HighFive::File &  h5file)
inherited

Definition at line 3736 of file miluphpc.cpp.

◆ distributionFromFile()

void Miluphpc::distributionFromFile ( const std::string &  filename)
inherited

Read initial/particle distribution file (in parallel)

Parameters
filenameinitial conditions file to be read

Definition at line 161 of file miluphpc.cpp.

◆ dynamicLoadBalancing()

void Miluphpc::dynamicLoadBalancing ( int  bins = 5000)
inherited

Pre-calculations for updateRangeApproximately.

Potential wrapper functions if more range determination functions come into exist.

Parameters
binsamount of bins the range will be subdivided

Definition at line 2885 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ energy()

real Miluphpc::energy ( )
privateinherited

Calculate the total amount of energy.

Returns
accumulated time of functions within

Definition at line 538 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ fixedLoadBalancing()

void Miluphpc::fixedLoadBalancing ( )
inherited

Load balancing via equidistant ranges.

Definition at line 2839 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ getMemoryInfo()

void Miluphpc::getMemoryInfo ( )
inherited

Definition at line 3760 of file miluphpc.cpp.

◆ gravity()

real Miluphpc::gravity ( )
inherited

Wrapper function for Gravity-related stuff.

The gravitational force computation corresponds to the Barnes-Hut method. First all relevant lowest domain list nodes belonging to another process are determined. Then the pseudo-particles and particles to be sent can be determined by checking the extended $\theta$-criterion for all nodes and all relevant lowest domain list nodes. This is done in an approach where particles are either marked to be sent or not in an extra array and afterwards copied to contiguous memory in order to send them via MPI. After receiving, the particles are inserted into the local tree, whereas first the pseudo-particles and then the particles are inserted. Since complete sub-trees are inserted, no new pseudo-particles are generated during this process. Afterwards, the actual force computation can be accomplished locally. Finally the inserted (pseudo-)particles are removed.

Returns
accumulated time for functions within

Definition at line 1241 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ integrate()

void Leapfrog::integrate ( int  step)
virtual

Implementation of the abstract integration method.

Parameters
stepIntegration step (number)

Implements Miluphpc.

Definition at line 21 of file leapfrog.cpp.

+ Here is the call graph for this function:

◆ numParticlesFromFile()

void Miluphpc::numParticlesFromFile ( const std::string &  filename)
inherited

Determine amount of particles (numParticles and numParticlesLocal) from initial file/particle distribution file.

Since the information of how many particles to be simulated is needed to allocate (the right amount of) memory, this function need to be called before distributionFromFile()

Parameters
filenameinitial conditions file to be read

Definition at line 123 of file miluphpc.cpp.

◆ parallel_gravity()

real Miluphpc::parallel_gravity ( )
privateinherited

Parallel version regarding computation of gravitational stuff.

Returns
accumulated time of functions within

Definition at line 1248 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ parallel_pseudoParticles()

real Miluphpc::parallel_pseudoParticles ( )
privateinherited

Parallel version regarding computation of pseudo-particles.

Returns
accumulated time of functions within

Definition at line 1061 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ parallel_sph()

real Miluphpc::parallel_sph ( )
privateinherited

Parallel version regarding computation of SPH-stuff.

Returns
accumulated time of functions within

Definition at line 2216 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ parallel_tree()

real Miluphpc::parallel_tree ( )
privateinherited

Parallel version regarding tree-stuff.

Returns
accumulated time of functions within

Definition at line 958 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ particles2file() [1/2]

real Miluphpc::particles2file ( const std::string &  filename,
int *  particleIndices,
int  length 
)
inherited

Write particles or rather particle's information/attributes to file. Particles to be outputted determined by particleIndices

Parameters
filenamefile name of output H5 file
particleIndicesparticle's indices to be outputted
lengthamount of particles to be outputted
Returns
accumulated time for functions within

Definition at line 3620 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ particles2file() [2/2]

real Miluphpc::particles2file ( int  step)
inherited

Write all of the information/particle distribution to H5 file

Parameters
stepsimulation step, used to name file
Returns
accumulated time for functions within

Definition at line 3409 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ prepareSimulation()

void Miluphpc::prepareSimulation ( )
inherited

Prepare the simulation, including.

  • loading the initial conditions
  • copying to device
  • computing the bounding boxes
  • initial ranges (load balancing)

Definition at line 235 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ pseudoParticles()

real Miluphpc::pseudoParticles ( )
inherited

Wrapper function for calculating pseudo-particles.

The pseudo-particle computation starts with determining the lowest domain list, this information is stored in a separate instance of the DomainList class continues with the local COM computation, the resetting of the domain list nodes that are not lowest domain list nodes, the preparation for the exchange by copying the information into contiguous memory, subsequent communication via MPI, the updating of the lowest domain list nodes and finally ends with the computation of the remaining domain list node pseudo-particles.

Returns
accumulated time for functions within

Definition at line 1054 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ removeParticles()

real Miluphpc::removeParticles ( )
inherited

Remove particles in dependence of some criterion.

Criterion can be specified in the config file. Removing particles can be accomplished by marking the particles to be removed via an extra array and using this one as key for sorting all the (pseudo-)particle arrays and finally deleting them by overwriting those entries with default values.

Returns
accumulated time of functions within

Definition at line 3040 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ reset()

real Miluphpc::reset ( )
privateinherited

Reset arrays, values, ...

This function embraces resetting the

  • pseudo-particles
  • child array
  • the boundaries
  • the domain list nodes
  • the neighbor list
  • and so forth ...

ensuring correct starting conditions.

Returns
accumulated time of functions within

Definition at line 586 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ rhs()

real Miluphpc::rhs ( int  step,
bool  selfGravity = true,
bool  assignParticlesToProcess = true 
)
inherited

Right Hand Side - all of the simulation without integration itself.

This method is used/called by the different integrators in order to integrate or execute the simulation. The force or acceleration computation and necessary steps to calculate them are encapsulated in this function whose individual components and interoperation are depicted in the following. This right-hand-side is called once or several times within a (sub-)integration step in dependence of the integration scheme. In addition, the function depends on whether the simulation runs single-node or multi-node, gravity and/or SPH are included into the simulation and some other conditions.

Parameters
stepsimulation step, regarding output
selfGravityapply self-gravity y/n (needed since gravity is decoupled)
assignParticlesToProcessassign particles to correct process y/n
Returns
accumulated time for functions within

Definition at line 294 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ sendParticles()

template<typename T >
integer Miluphpc::sendParticles ( T *  sendBuffer,
T *  receiveBuffer,
integer *  sendLengths,
integer *  receiveLengths 
)
privateinherited

Send particles/Exchange particles among MPI processes.

Template Parameters
T
Parameters
sendBuffer
receiveBuffer
sendLengths
receiveLengths
Returns

Definition at line 3174 of file miluphpc.cpp.

◆ sendParticlesEntry()

template<typename T >
integer Miluphpc::sendParticlesEntry ( integer *  sendLengths,
integer *  receiveLengths,
T *  entry,
T *  entryBuffer,
T *  copyBuffer 
)
privateinherited

Send particles/Exchange particles among MPI processes.

Template Parameters
T
Parameters
sendLengths
receiveLengths
entry
entryBuffer
copyBuffer
Returns

Definition at line 3327 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ sph()

real Miluphpc::sph ( )
inherited

Wrapper function for SPH-related stuff

Similar to the gravitational force computation is the SPH part. The relevant lowest domain list nodes are determined in order to decide which particles need to be sent, which is also done in an analogous approach in the gravitational part. However, in this case only particles are sent and after exchanging those particles, the insertion includes generation of new pseudo-particles similar to building or rather extending the local tree. With inserted particles from the other processes the FRNN search can be done and subsequent computation of the density, speed of sound and pressure realized. The corresponding particle properties are sent again, so that the actual force computation can be performed. Finally, the inserted and generated (pseudo-)particles are removed.

Returns
accumulated time for functions within

Definition at line 2207 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ tree()

real Miluphpc::tree ( )
inherited

Wrapper function for building the tree (and domain tree).

To build the parallel tree it is necessary to determine the domain list or common coarse tree nodes which are derived from the ranges. After that, the tree is built locally and finally all the domain list nodes or rather those attributes as indices are either saved to an instance of the omainList class or if they do not exist after building the tree created and afterwards assigned.

Returns
accumulated time for functions within

Definition at line 951 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ updateRange()

void Miluphpc::updateRange ( int  aimedParticlesPerProcess)
inherited

Update the range in dependence on number of (MPI) processes and aimed particles per process.

Counting the particles per process, calculating the aimed particles per process as quotient of total number of particles and number of processes, determining the particle keys, sorting them and determining the new ranges as the indices of the multiplies of the aimed number of particles per process.

Parameters
aimedParticlesPerProcessoptimal number of particles per process

Definition at line 2912 of file miluphpc.cpp.

+ Here is the call graph for this function:

◆ updateRangeApproximately()

void Miluphpc::updateRangeApproximately ( int  aimedParticlesPerProcess,
int  bins = 5000 
)
inherited

Update the ranges (approximately and dynamically).

A histogram is generated with the amount of bins given by bins. The possible range is distributed accordingly and particles are sorted in dependence of their key into, in order to extract the ranges.

Parameters
aimedParticlesPerProcessaimed amount particles per process
binsamount of bins the range will be subdivided

Definition at line 2988 of file miluphpc.cpp.

+ Here is the call graph for this function:

Member Data Documentation

◆ buffer

HelperHandler* Miluphpc::buffer
inherited

buffer instance

Todo:
revise buffer handling buffer instance

Definition at line 325 of file miluphpc.h.

◆ curveType

Curve::Type Miluphpc::curveType
inherited

Space-filling curve type to be used (Lebesgue or Hilbert)

Definition at line 277 of file miluphpc.h.

◆ d_mutex

integer* Miluphpc::d_mutex
inherited

Definition at line 307 of file miluphpc.h.

◆ domainListHandler

DomainListHandler* Miluphpc::domainListHandler
inherited

Instance to handle the DomainList instance on device and host.

Definition at line 315 of file miluphpc.h.

◆ h_searchRadius

real Miluphpc::h_searchRadius
inherited

search radius for SPH (MPI-process overarching) neighbor search

Definition at line 215 of file miluphpc.h.

◆ integratedParticles

IntegratedParticleHandler* Miluphpc::integratedParticles
inherited

Instance(s) to handle the IntegratedParticles instance(s) on device and host

Memory for this may or may not be allocated within child classes. This depends on whether instances are needed for cache particle information e.g. for predictor-corrector integrators.

Definition at line 305 of file miluphpc.h.

◆ kernelHandler

SPH::KernelHandler Miluphpc::kernelHandler
inherited

Instance to handle the SPH Kernel instance on device and host.

Definition at line 280 of file miluphpc.h.

◆ lowestDomainListHandler

DomainListHandler* Miluphpc::lowestDomainListHandler
inherited

Instance to handle the (lowest) DomainList instance on device and host.

Definition at line 317 of file miluphpc.h.

◆ materialHandler

MaterialHandler* Miluphpc::materialHandler
inherited

Instance to handle Materials instances on device and host.

Definition at line 319 of file miluphpc.h.

◆ numNodes

integer Miluphpc::numNodes
inherited

number of nodes (to be allocated)

Definition at line 296 of file miluphpc.h.

◆ numParticles

integer Miluphpc::numParticles
inherited

number of particles (to be allocated)

Definition at line 285 of file miluphpc.h.

◆ numParticlesLocal

integer Miluphpc::numParticlesLocal
inherited

number of particles currently living on the (MPI) process

Temporarily are more particles on the process, e.g. for gravitational and SPH (force) calculations.

Definition at line 294 of file miluphpc.h.

◆ particleHandler

ParticleHandler* Miluphpc::particleHandler
inherited

Instance to handle the Particles instance on device and host.

Definition at line 309 of file miluphpc.h.

◆ profiler

H5Profiler& Miluphpc::profiler = H5Profiler::getInstance("log/performance.h5")
inherited

H5 profiler instance.

Definition at line 274 of file miluphpc.h.

◆ simulationParameters

SimulationParameters Miluphpc::simulationParameters
inherited

buffer (need for revising)

buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) buffer (need for revising) collected information required to set up the simulation

Definition at line 351 of file miluphpc.h.

◆ simulationTimeHandler

SimulationTimeHandler* Miluphpc::simulationTimeHandler
inherited

Instance to handle the SimulationTime instances on device and host.

Definition at line 282 of file miluphpc.h.

◆ subDomainKeyTreeHandler

SubDomainKeyTreeHandler* Miluphpc::subDomainKeyTreeHandler
inherited

Instance to handle the SubDomainKeyTree instance on device and host.

Definition at line 311 of file miluphpc.h.

◆ subStep

int Miluphpc::subStep
inherited

current sub-step (there are possibly more sub-steps within a step!)

Definition at line 212 of file miluphpc.h.

◆ sumParticles

integer Miluphpc::sumParticles
inherited

(real) number of particles on all processes

Definition at line 287 of file miluphpc.h.

◆ totalAngularMomentum

real Miluphpc::totalAngularMomentum
inherited

total angular momentum

Definition at line 221 of file miluphpc.h.

◆ totalEnergy

real Miluphpc::totalEnergy
inherited

total energy

Definition at line 218 of file miluphpc.h.

◆ treeHandler

TreeHandler* Miluphpc::treeHandler
inherited

Instance to handle the Tree instance on device and host.

Definition at line 313 of file miluphpc.h.


The documentation for this class was generated from the following files:
  • include/integrator/leapfrog.h
  • src/integrator/leapfrog.cpp

milupHPC - Leapfrog Class Reference
Generated on Wed Aug 31 2022 12:16:53 by Doxygen 1.9.3