Gromacs  2016.5
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
pbc-simd.h File Reference
#include "config.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/simd/simd.h"
+ Include dependency graph for pbc-simd.h:
+ This graph shows which files directly or indirectly include this file:

Description

This file contains a definition, declaration and inline function for SIMD accelerated PBC calculations.

Author
Berk Hess hess@.nosp@m.kth..nosp@m.se

Functions

void set_pbc_simd (const t_pbc *pbc, real *pbc_simd)
 Set the SIMD PBC data from a normal t_pbc struct. More...
 
static void gmx_simdcall pbc_correct_dx_simd (SimdReal *dx, SimdReal *dy, SimdReal *dz, const real *pbc_simd)
 Correct SIMD distance vector *dx,*dy,*dz for PBC using SIMD. More...
 
static void gmx_simdcall pbc_dx_aiuc (const real *pbc_simd, const SimdReal *x1, const SimdReal *x2, SimdReal *dx)
 Calculates the PBC corrected distance between SIMD coordinates. More...
 

Function Documentation

static void gmx_simdcall pbc_correct_dx_simd ( SimdReal dx,
SimdReal dy,
SimdReal dz,
const real pbc_simd 
)
inlinestatic

Correct SIMD distance vector *dx,*dy,*dz for PBC using SIMD.

For rectangular boxes all returned distance vectors are the shortest. For triclinic boxes only distances up to half the smallest box diagonal element are guaranteed to be the shortest. This means that distances from 0.5/sqrt(2) times a box vector length (e.g. for a rhombic dodecahedron) can use a more distant periodic image. Note that this routine always does PBC arithmetic, even for dimensions without PBC. But on modern processors the overhead of this, often called, routine should be low. On e.g. Intel Haswell/Broadwell it takes 8 cycles.

static void gmx_simdcall pbc_dx_aiuc ( const real pbc_simd,
const SimdReal x1,
const SimdReal x2,
SimdReal dx 
)
inlinestatic

Calculates the PBC corrected distance between SIMD coordinates.

Parameters
pbc_simdSIMD formatted PBC information
x1Packed coordinates of atom1, size 3*GMX_SIMD_REAL_WIDTH
x2Packed coordinates of atom2, size 3*GMX_SIMD_REAL_WIDTH
dxThe PBC corrected distance x1 - x2

This routine only returns the shortest distance correctd for PBC when all atoms are in the unit-cell (aiuc). This is the SIMD equivalent of the scalar version declared in pbc.h.

void set_pbc_simd ( const t_pbc pbc,
real pbc_simd 
)

Set the SIMD PBC data from a normal t_pbc struct.

Parameters
pbcType of periodic boundary, NULL can be passed for then no PBC will be used.
pbc_simdPointer to aligned memory with (DIM + DIM*(DIM+1)/2) GMX_SIMD_REAL_WIDTH reals describing the box vectors unrolled by GMX_SIMD_REAL_WIDTH. These are sorted in a slightly non-standard order so that we always issue the memory loads in order (to improve prefetching) in pbc_correct_dx_simd(). The order is inv_bzz, bzx, bzy, bzz, inv_byy, byx, byy, inv_bxx, and bxx.