Computation-reordering and cache friendly computation
In this example we show how reordering the data can significantly improve the computation speed. In order to do this we will re-work the molecular dynamic example.
Calculate forces
This function is the same as the molecular dynamic example with few changes:
- The function now take as argument CellList_hilb instead of CellList
template<
typename CellList>
void calc_forces(
vector_dist<3,
double,
aggregate<
double[3],
double[3]> > & vd,
CellList & NN,
double sigma12,
double sigma6)
{
Class for FAST cell list implementation.
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
- We get an iterator from the Cell list instead that from the vector
auto it2 = NN.getIterator();
- See also
- Calculate forces
Calculate energy
In this function is the same as the molecular dynamic example with few changes:
- The function now take as argument CellList_hilb instead of CellList
template<
typename CellList>
void calc_forces(
vector_dist<3,
double,
aggregate<
double[3],
double[3]> > & vd,
CellList & NN,
double sigma12,
double sigma6)
{
- We get an iterator from the Cell list instead that from the vector
auto it2 = NN.getIterator();
The difference in doing this is that now we iterate on particles in a smarter way. We will explain more in detail later in the example
- See also
- Calculate energy
Initialization
The initialization is the same as the molecular dynamic example. The differences are in the parameters. We will use a bigger system, with more particles. The delta time for integration is chosen in order to keep the system stable.
- See also
- Initialization
double dt = 0.0001;
double r_cut = 0.03;
double sigma = r_cut/3.0;
double sigma12 = pow(sigma,12);
double sigma6 = pow(sigma,6);
openfpm_init(&argc,&argv);
size_t sz[3] = {40,40,40};
size_t bc[3]={PERIODIC,PERIODIC,PERIODIC};
This class represent an N-dimensional box.
Implementation of VCluster class.
Implementation of 1-D std::vector like structure.
Particles on a grid like position
Here we place the particles on a grid like manner
- See also
- Particles on a grid like position
auto it = vd.getGridIterator(sz);
while (it.isNext())
{
vd.add();
auto key = it.get();
vd.getLastPos()[0] = key.get(0) * it.getSpacing(0);
vd.getLastPos()[1] = key.get(1) * it.getSpacing(1);
vd.getLastPos()[2] = key.get(2) * it.getSpacing(2);
vd.template getLastProp<velocity>()[0] = 0.0;
vd.template getLastProp<velocity>()[1] = 0.0;
vd.template getLastProp<velocity>()[2] = 0.0;
vd.template getLastProp<force>()[0] = 0.0;
vd.template getLastProp<force>()[1] = 0.0;
vd.template getLastProp<force>()[2] = 0.0;
++it;
}
Molecular dynamic steps
Here we do 30000 MD steps using verlet integrator the cycle is the same as the molecular dynamic example. with the following changes.
Cell lists
Instead of getting the normal cell list we get an hilbert curve cell-list. Such cell list has a function called getIterator used inside the function calc_forces and calc_energy that iterate across all the particles but in a smart-way. In practice given an r-cut a cell-list is constructed with the provided spacing. Suppose to have a cell-list \( m \times n \), an hilbert curve \( 2^k \times 2^k \) is contructed with \( k = ceil(log_2(max(m,n))) \). Cell-lists are explored according to this Hilbert curve, If a cell does not exist is simply skipped.
+------+------+------+------+ Example of Hilbert curve running on a 3 x 3 Cell
| | | | | An hilbert curve of k = ceil(log_2(3)) = 4
| X+---->X | X +---> X |
| ^ | + | ^ | + |
***|******|******|****---|--+ *******
* + | v | + * v | * *
* 7 | 8+---->9 * X | * * = Domain
* ^ | | * + | * *
*--|-----------------*---|--+ *******
* + | | * v |
* 4<----+5 | 6<---+ X |
* | ^ | + * |
*---------|-------|--*------+
* | + | v * |
* 1+---->2 | 3+---> X |
* | | * |
**********************------+
this mean that we will iterate the following cells
1,2,5,4,7,8,9,6,3
Suppose now that the particles are ordered like described
Particles id Cell
0 1
1 7
2 8
3 1
4 9
5 9
6 6
7 7
8 3
9 2
10 4
11 3
The iterator of the cell-list will explore the particles in the following way
Cell 1 2 5 4 7 8 9 6 3
| | | | | | | | | |
0,3,9,,10,1,7,2,4,5,6,8
*
We cannot explain here what is a cache, but in practice is a fast memory in the CPU able to store chunks of memory. The cache in general is much smaller than RAM, but the big advantage is its speed. Retrieve data from the cache is much faster than RAM. Unfortunately the factors that determine what is on cache and what is not are multiples: Type of cache, algorithm ... . Qualitatively all caches will tend to load chunks of data that you read multiple-time, or chunks of data that probably you will read based on pattern analysis. A small example is a linear memory copy where you read consecutively memory and you write on consecutive memory. Modern CPU recognize such pattern and decide to load on cache the consecutive memory before you actually require it.
Iterating the vector in the way described above has the advantage that when we do computation on particles and its neighborhood with the sequence described above it will happen that:
- If to process a particle A we read some neighborhood particles to process the next particle A+1 we will probably read most of the previous particles.
In order to show in practice what happen we first show the graph when we do not reorder
The measure has oscillation but we see an asymptotic behavior from 0.04 in the initial condition to 0.124 . Below we show what happen when we use iterator from the Cell list hilbert
In cases where particles does not move or move very slowly consider to use data-reordering, because it can give 8-10% speedup
- See also
- Data reordering, computation-reordering and cache friendly computation
Timers
In order to collect the time of the force calculation we insert two timers around the function calc_force. The overall performance is instead calculated with another timer around the time stepping
if (i % 10 == 0)
Class for cpu time benchmarking.
void start()
Start the timer.
if (i % 10 == 0)
{
x.add(i);
}
void stop()
Stop the timer.
double getwct()
Return the elapsed real time.
- See also
- Molecular dynamic steps
Plotting graphs
After we terminate the MD steps our vector x contains at which iteration we benchmark the force calculation time, while y contains the measured time at that time-step. We can produce a graph X Y
- Note
- The graph produced is an svg graph that can be view with a browser. From the browser we can also easily save the graph into pure svg format
options.
title = std::string(
"Force calculation time");
options.
yAxis = std::string(
"Time");
options.
xAxis = std::string(
"iteration");
cg.
write(
"gc_plot2_out.html");
Small class to produce graph with Google chart in HTML.
void write(std::string file)
It write the graphs on file in html format using Google charts.
void AddLinesGraph(openfpm::vector< X > &x, openfpm::vector< Y > &y, const GCoptions &opt)
Add a simple lines graph.
std::string xAxis
X axis name.
size_t lineWidth
Width of the line.
std::string title
Title of the chart.
std::string yAxis
Y axis name.
Finalize
At the very end of the program we have always to de-initialize the library
Full code
#include "Vector/vector_dist.hpp"
#include "Decomposition/CartDecomposition.hpp"
#include "data_type/aggregate.hpp"
#include "Plot/GoogleChart.hpp"
#include "Plot/util.hpp"
#include "timer.hpp"
constexpr int velocity = 0;
constexpr int force = 1;
template<
typename CellList>
void calc_forces(
vector_dist<3,
double,
aggregate<
double[3],
double[3]> > & vd,
CellList & NN,
double sigma12,
double sigma6)
{
vd.updateCellList(NN);
auto it2 = NN.getIterator();
while (it2.isNext())
{
vd.template getProp<force>(p)[0] = 0.0;
vd.template getProp<force>(p)[1] = 0.0;
vd.template getProp<force>(p)[2] = 0.0;
auto Np = NN.template getNNIterator<NO_CHECK>(NN.getCell(vd.getPos(p)));
while (Np.isNext())
{
auto q = Np.get();
if (q == p) {++Np; continue;};
double rn = norm2(r);
Point<3,double> f = 24.0*(2.0 *sigma12 / (rn*rn*rn*rn*rn*rn*rn) - sigma6 / (rn*rn*rn*rn)) * r;
vd.template getProp<force>(p)[0] += f.get(0);
vd.template getProp<force>(p)[1] += f.get(1);
vd.template getProp<force>(p)[2] += f.get(2);
++Np;
}
++it2;
}
}
template<
typename CellList>
double calc_energy(
vector_dist<3,
double,
aggregate<
double[3],
double[3]> > & vd,
CellList & NN,
double sigma12,
double sigma6)
{
double E = 0.0;
vd.updateCellList(NN);
auto it2 = NN.getIterator();
while (it2.isNext())
{
vd.template getProp<force>(p)[0] = 0.0;
vd.template getProp<force>(p)[1] = 0.0;
vd.template getProp<force>(p)[2] = 0.0;
auto Np = NN.template getNNIterator<NO_CHECK>(NN.getCell(vd.getPos(p)));
while (Np.isNext())
{
auto q = Np.get();
if (q == p) {++Np; continue;};
double rn = norm2(xp - xq);
E += 4.0 * ( sigma12 / (rn*rn*rn*rn*rn*rn) - sigma6 / ( rn*rn*rn) );
++Np;
}
E += (vd.template getProp<velocity>(p)[0]*vd.template getProp<velocity>(p)[0] +
vd.template getProp<velocity>(p)[1]*vd.template getProp<velocity>(p)[1] +
vd.template getProp<velocity>(p)[2]*vd.template getProp<velocity>(p)[2]) / 2;
++it2;
}
return E;
}
int main(int argc, char* argv[])
{
double dt = 0.0001;
double r_cut = 0.03;
double sigma = r_cut/3.0;
double sigma12 = pow(sigma,12);
double sigma6 = pow(sigma,6);
openfpm_init(&argc,&argv);
size_t sz[3] = {40,40,40};
size_t bc[3]={PERIODIC,PERIODIC,PERIODIC};
auto it = vd.getGridIterator(sz);
while (it.isNext())
{
vd.add();
auto key = it.get();
vd.getLastPos()[0] = key.get(0) * it.getSpacing(0);
vd.getLastPos()[1] = key.get(1) * it.getSpacing(1);
vd.getLastPos()[2] = key.get(2) * it.getSpacing(2);
vd.template getLastProp<velocity>()[0] = 0.0;
vd.template getLastProp<velocity>()[1] = 0.0;
vd.template getLastProp<velocity>()[2] = 0.0;
vd.template getLastProp<force>()[0] = 0.0;
vd.template getLastProp<force>()[1] = 0.0;
vd.template getLastProp<force>()[2] = 0.0;
++it;
}
auto NN = vd.getCellList_hilb(r_cut);
calc_forces(vd,NN,sigma12,sigma6);
unsigned long int f = 0;
#ifndef TEST_RUN
size_t Nstep = 30000;
#else
size_t Nstep = 300;
#endif
for (size_t i = 0; i < Nstep ; i++)
{
auto it3 = vd.getDomainIterator();
while (it3.isNext())
{
auto p = it3.get();
vd.template getProp<velocity>(p)[0] += 0.5*dt*vd.template getProp<force>(p)[0];
vd.template getProp<velocity>(p)[1] += 0.5*dt*vd.template getProp<force>(p)[1];
vd.template getProp<velocity>(p)[2] += 0.5*dt*vd.template getProp<force>(p)[2];
vd.getPos(p)[0] += vd.template getProp<velocity>(p)[0]*dt;
vd.getPos(p)[1] += vd.template getProp<velocity>(p)[1]*dt;
vd.getPos(p)[2] += vd.template getProp<velocity>(p)[2]*dt;
++it3;
}
vd.map();
vd.template ghost_get<>();
if (i % 10 == 0)
calc_forces(vd,NN,sigma12,sigma6);
if (i % 10 == 0)
{
x.add(i);
}
auto it4 = vd.getDomainIterator();
while (it4.isNext())
{
auto p = it4.get();
vd.template getProp<velocity>(p)[0] += 0.5*dt*vd.template getProp<force>(p)[0];
vd.template getProp<velocity>(p)[1] += 0.5*dt*vd.template getProp<force>(p)[1];
vd.template getProp<velocity>(p)[2] += 0.5*dt*vd.template getProp<force>(p)[2];
++it4;
}
if (i % 100 == 0)
{
vd.deleteGhost();
vd.write("particles_",f);
vd.ghost_get<>();
double energy = calc_energy(vd,NN,sigma12,sigma6);
auto & vcl = create_vcluster();
vcl.sum(energy);
vcl.execute();
if (vcl.getProcessUnitID() == 0)
std::cout << std::endl << "Energy: " << energy << std::endl;
f++;
}
}
std::cout <<
"Performance: " << time2.
getwct() << std::endl;
options.
title = std::string(
"Force calculation time");
options.
yAxis = std::string(
"Time");
options.
xAxis = std::string(
"iteration");
cg.
write(
"gc_plot2_out.html");
openfpm_finalize();
}
auto get(size_t cell, size_t ele) -> decltype(this->Mem_type::get(cell, ele))
Get an element in the cell.
This class implement the point shape in an N-dimensional space.