#include "leapfrog.h"
|
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...
|
|
Definition at line 16 of file leapfrog.h.
◆ Leapfrog()
Constructor.
- Parameters
-
simulationParameters | Simulation parameters/settings. |
Definition at line 9 of file leapfrog.cpp.
◆ ~Leapfrog()
◆ afterIntegrationStep()
void Miluphpc::afterIntegrationStep |
( |
| ) |
|
|
inherited |
Function to be called after successful integration (step)
Definition at line 471 of file miluphpc.cpp.
◆ 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.
◆ 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
-
U | type of the arrays determining the sorting |
T | type of the arrays to be sorted |
- Parameters
-
sortArray | array to be sorted (rather sorting behaviour should be determined) |
sortedArray | sorted array (result of sorting) |
entry | (relevant) array to be sorted |
temp | buffer needed for sorting process |
- Returns
- accumulated time of functions within
Definition at line 944 of file miluphpc.cpp.
◆ 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.
◆ 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.
◆ configParameters2file()
real Miluphpc::configParameters2file |
( |
HighFive::File & |
h5file | ) |
|
|
inherited |
◆ distributionFromFile()
void Miluphpc::distributionFromFile |
( |
const std::string & |
filename | ) |
|
|
inherited |
Read initial/particle distribution file (in parallel)
- Parameters
-
filename | initial 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
-
bins | amount of bins the range will be subdivided |
Definition at line 2885 of file miluphpc.cpp.
◆ 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.
◆ fixedLoadBalancing()
void Miluphpc::fixedLoadBalancing |
( |
| ) |
|
|
inherited |
Load balancing via equidistant ranges.
Definition at line 2839 of file miluphpc.cpp.
◆ getMemoryInfo()
void Miluphpc::getMemoryInfo |
( |
| ) |
|
|
inherited |
◆ 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.
◆ integrate()
void Leapfrog::integrate |
( |
int |
step | ) |
|
|
virtual |
Implementation of the abstract integration method.
- Parameters
-
step | Integration step (number) |
Implements Miluphpc.
Definition at line 21 of file leapfrog.cpp.
◆ 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
-
filename | initial 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.
◆ 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.
◆ 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.
◆ 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.
◆ 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
-
filename | file name of output H5 file |
particleIndices | particle's indices to be outputted |
length | amount of particles to be outputted |
- Returns
- accumulated time for functions within
Definition at line 3620 of file miluphpc.cpp.
◆ particles2file() [2/2]
real Miluphpc::particles2file |
( |
int |
step | ) |
|
|
inherited |
Write all of the information/particle distribution to H5 file
- Parameters
-
step | simulation step, used to name file |
- Returns
- accumulated time for functions within
Definition at line 3409 of file miluphpc.cpp.
◆ 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.
◆ 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.
◆ 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.
◆ reset()
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.
◆ 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
-
step | simulation step, regarding output |
selfGravity | apply self-gravity y/n (needed since gravity is decoupled) |
assignParticlesToProcess | assign particles to correct process y/n |
- Returns
- accumulated time for functions within
Definition at line 294 of file miluphpc.cpp.
◆ 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
-
- 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
-
- Parameters
-
sendLengths | |
receiveLengths | |
entry | |
entryBuffer | |
copyBuffer | |
- Returns
Definition at line 3327 of file miluphpc.cpp.
◆ sph()
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.
◆ tree()
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.
◆ 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
-
aimedParticlesPerProcess | optimal number of particles per process |
Definition at line 2912 of file miluphpc.cpp.
◆ 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
-
aimedParticlesPerProcess | aimed amount particles per process |
bins | amount of bins the range will be subdivided |
Definition at line 2988 of file miluphpc.cpp.
◆ buffer
buffer instance
- Todo:
- revise buffer handling buffer instance
Definition at line 325 of file miluphpc.h.
◆ curveType
Space-filling curve type to be used (Lebesgue or Hilbert)
Definition at line 277 of file miluphpc.h.
◆ d_mutex
◆ domainListHandler
◆ 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
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
◆ lowestDomainListHandler
◆ materialHandler
Instance to handle Materials
instances on device and host.
Definition at line 319 of file miluphpc.h.
◆ numNodes
number of nodes (to be allocated)
Definition at line 296 of file miluphpc.h.
◆ numParticles
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
◆ profiler
◆ simulationParameters
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
◆ subDomainKeyTreeHandler
◆ subStep
current sub-step (there are possibly more sub-steps within a step!)
Definition at line 212 of file miluphpc.h.
◆ sumParticles
(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 |
◆ treeHandler
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: