Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Macros | Enumerations | Functions | Variables
pbc.cpp File Reference
#include "gmxpre.h"
#include "pbc.h"
#include <cmath>
#include <algorithm>
#include "gromacs/math/functions.h"
#include "gromacs/math/units.h"
#include "gromacs/math/utilities.h"
#include "gromacs/math/vec.h"
#include "gromacs/math/vecdump.h"
#include "gromacs/pbcutil/ishift.h"
#include "gromacs/pbcutil/mshift.h"
#include "gromacs/topology/topology.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/smalloc.h"
+ Include dependency graph for pbc.cpp:

Description

Implements routines in pbc.h.

Utility functions for handling periodic boundary conditions. Mainly used in analysis tools.

Macros

#define BOX_MARGIN   1.0010
 Margin factor for error message.
 
#define BOX_MARGIN_CORRECT   1.0005
 Margin correction if the box is too skewed.
 

Enumerations

enum  {
  epbcdxRECTANGULAR = 1, epbcdxTRICLINIC, epbcdx2D_RECT, epbcdx2D_TRIC,
  epbcdx1D_RECT, epbcdx1D_TRIC, epbcdxSCREW_RECT, epbcdxSCREW_TRIC,
  epbcdxNOPBC, epbcdxUNSUPPORTED
}
 

Functions

int ePBC2npbcdim (int ePBC)
 Returns the number of dimensions that use pbc. More...
 
void dump_pbc (FILE *fp, t_pbc *pbc)
 Dump the contents of the pbc structure to the file. More...
 
const char * check_box (int ePBC, const matrix box)
 Check the box for consistency. More...
 
void matrix_convert (matrix box, const rvec vec, const rvec angleInDegrees)
 Creates box matrix from edge lengths and angles. More...
 
real max_cutoff2 (int ePBC, const matrix box)
 Compute the maximum cutoff for the box. More...
 
int guess_ePBC (const matrix box)
 Guess PBC typr. More...
 
static int correct_box_elem (FILE *fplog, int step, tensor box, int v, int d)
 Check if the box still obeys the restrictions, if not, correct it.
 
gmx_bool correct_box (FILE *fplog, int step, tensor box, t_graph *graph)
 Corrects the box if necessary. More...
 
static void low_set_pbc (t_pbc *pbc, int ePBC, const ivec dd_pbc, const matrix box)
 Do the real arithmetic for filling the pbc struct.
 
void set_pbc (t_pbc *pbc, int ePBC, const matrix box)
 Initiate the periodic boundary condition algorithms. More...
 
t_pbcset_pbc_dd (t_pbc *pbc, int ePBC, const ivec domdecCells, gmx_bool bSingleDir, const matrix box)
 Initiate the periodic boundary condition algorithms. More...
 
void pbc_dx (const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
 Compute distance with PBC. More...
 
int pbc_dx_aiuc (const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
 Compute distance vector for simple PBC types. More...
 
void pbc_dx_d (const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx)
 Compute distance vector in double precision. More...
 
void calc_shifts (const matrix box, rvec shift_vec[])
 Computes shift vectors. More...
 
void calc_box_center (int ecenter, const matrix box, rvec box_center)
 Calculates the center of the box. More...
 
void calc_triclinic_images (const matrix box, rvec img[])
 Calculates the NTRICIMG box images. More...
 
void calc_compact_unitcell_vertices (int ecenter, const matrix box, rvec vert[])
 Calculates the NCUCVERT vertices of a compact unitcell. More...
 
int * compact_unitcell_edges ()
 Compute unitcell edges. More...
 
void put_atoms_in_box (int ePBC, const matrix box, gmx::ArrayRef< gmx::RVec > x)
 Put atoms inside the simulations box. More...
 
void put_atoms_in_box_omp (int ePBC, const matrix box, gmx::ArrayRef< gmx::RVec > x, int nth)
 Parallellizes put_atoms_in_box() More...
 
void put_atoms_in_triclinic_unitcell (int ecenter, const matrix box, gmx::ArrayRef< gmx::RVec > x)
 Put atoms inside triclinic box. More...
 
void put_atoms_in_compact_unitcell (int ePBC, int ecenter, const matrix box, gmx::ArrayRef< gmx::RVec > x)
 Put atoms inside the unitcell. More...
 
static void low_do_pbc_mtop (FILE *fplog, int ePBC, const matrix box, const gmx_mtop_t *mtop, rvec x[], gmx_bool bFirst)
 Make molecules whole by shifting positions. More...
 
void do_pbc_first_mtop (FILE *fplog, int ePBC, const matrix box, const gmx_mtop_t *mtop, rvec x[])
 Make all molecules whole by shifting positions. More...
 
void do_pbc_mtop (int ePBC, const matrix box, const gmx_mtop_t *mtop, rvec x[])
 Make molecules consisting of multiple charge groups whole by shifting positions. More...
 

Variables

const char * epbc_names [epbcNR+1] = { "xyz", "no", "xy", "screw", nullptr }
 Strings corresponding to epbc enum values.
 
static gmx_bool bWarnedGuess = false
 Set to true if warning has been printed.
 

Function Documentation

void calc_box_center ( int  ecenter,
const matrix  box,
rvec  box_center 
)

Calculates the center of the box.

See the description for the enum ecenter above.

Parameters
[in]ecenterDescription of center type
[in]boxThe simulation box
[out]box_centerThe center of the box
void calc_compact_unitcell_vertices ( int  ecenter,
const matrix  box,
rvec  vert[] 
)

Calculates the NCUCVERT vertices of a compact unitcell.

Parameters
[in]ecenterThe center type
[in]boxThe simulation box
[out]vertThe vertices
void calc_shifts ( const matrix  box,
rvec  shift_vec[] 
)

Computes shift vectors.

This routine calculates ths shift vectors necessary to use the neighbor searching routine.

Parameters
[in]boxThe simulation box
[out]shift_vecThe shifting vectors
void calc_triclinic_images ( const matrix  box,
rvec  img[] 
)

Calculates the NTRICIMG box images.

Parameters
[in]boxThe simulation box
[out]imgThe triclinic box images
const char* check_box ( int  ePBC,
const matrix  box 
)

Check the box for consistency.

Parameters
[in]ePBCThe pbc identifier
[in]boxThe box matrix
Returns
NULL if the box is supported by Gromacs. Otherwise returns a string with the problem. When ePBC=-1, the type of pbc is guessed from the box matrix.
int* compact_unitcell_edges ( )

Compute unitcell edges.

Returns
an array of unitcell edges of length NCUCEDGE*2, this is an index in vert[], which is calculated by calc_unitcell_vertices. The index consists of NCUCEDGE pairs of vertex indices. The index does not change, so it needs to be retrieved only once.
gmx_bool correct_box ( FILE *  fplog,
int  step,
tensor  box,
struct t_graph *  graph 
)

Corrects the box if necessary.

Checks for un-allowed box angles and corrects the box and the integer shift vectors in the graph (if graph!=NULL) if necessary.

Parameters
[in]fplogFile for debug output
[in]stepThe MD step number
[in]boxThe simulation cell
[in]graphInformation about molecular connectivity
Returns
TRUE when the box was corrected.
void do_pbc_first_mtop ( FILE *  fplog,
int  ePBC,
const matrix  box,
const gmx_mtop_t *  mtop,
rvec  x[] 
)

Make all molecules whole by shifting positions.

Parameters
[in]fplogLog file
[in]ePBCThe PBC type
[in]boxThe simulation box
[in]mtopSystem topology definition
[in,out]xThe coordinates of the atoms
void do_pbc_mtop ( int  ePBC,
const matrix  box,
const gmx_mtop_t *  mtop,
rvec  x[] 
)

Make molecules consisting of multiple charge groups whole by shifting positions.

Parameters
[in]ePBCThe PBC type
[in]boxThe simulation box
[in]mtopSystem topology definition
[in,out]xThe coordinates of the atoms
void dump_pbc ( FILE *  fp,
t_pbc pbc 
)

Dump the contents of the pbc structure to the file.

Parameters
[in]fpThe file pointer to write to
[in]pbcThe periodic boundary condition information structure
int ePBC2npbcdim ( int  ePBC)

Returns the number of dimensions that use pbc.

Parameters
[in]ePBCThe periodic boundary condition type
Returns
the number of dimensions that use pbc, starting at X
int guess_ePBC ( const matrix  box)

Guess PBC typr.

Guesses the type of periodic boundary conditions using the box

Parameters
[in]boxThe box matrix
Returns
The pbc identifier
static void low_do_pbc_mtop ( FILE *  fplog,
int  ePBC,
const matrix  box,
const gmx_mtop_t *  mtop,
rvec  x[],
gmx_bool  bFirst 
)
static

Make molecules whole by shifting positions.

Parameters
[in]fplogLog file
[in]ePBCThe PBC type
[in]boxThe simulation box
[in]mtopSystem topology definition
[in,out]xThe coordinates of the atoms
[in]bFirstSpecifier for first-time PBC removal
void matrix_convert ( matrix  box,
const rvec  vec,
const rvec  angleInDegrees 
)

Creates box matrix from edge lengths and angles.

Parameters
[in,out]boxThe box matrix
[in]vecThe edge lengths
[in]angleInDegreesThe angles
real max_cutoff2 ( int  ePBC,
const matrix  box 
)

Compute the maximum cutoff for the box.

Returns the square of the maximum cut-off allowed for the box, taking into account that the grid neighborsearch code and pbc_dx only check combinations of single box-vector shifts.

Parameters
[in]ePBCThe pbc identifier
[in]boxThe box matrix
Returns
the maximum cut-off.
void pbc_dx ( const t_pbc pbc,
const rvec  x1,
const rvec  x2,
rvec  dx 
)

Compute distance with PBC.

Calculate the correct distance vector from x2 to x1 and put it in dx. set_pbc must be called before ever calling this routine.

Note that for triclinic boxes that do not obey the GROMACS unit-cell restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.

Parameters
[in,out]pbcThe pbc information structure
[in]x1Coordinates for particle 1
[in]x2Coordinates for particle 2
[out]dxDistance vector
int pbc_dx_aiuc ( const t_pbc pbc,
const rvec  x1,
const rvec  x2,
rvec  dx 
)

Compute distance vector for simple PBC types.

Calculate the correct distance vector from x2 to x1 and put it in dx, This function can only be used when all atoms are in the rectangular or triclinic unit-cell. set_pbc_dd or set_pbc must be called before ever calling this routine.

Parameters
[in,out]pbcThe pbc information structure
[in]x1Coordinates for particle 1
[in]x2Coordinates for particle 2
[out]dxDistance vector
Returns
the ishift required to shift x1 at closest distance to x2; i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx (see calc_shifts below on how to obtain shift_vec)
void pbc_dx_d ( const t_pbc pbc,
const dvec  x1,
const dvec  x2,
dvec  dx 
)

Compute distance vector in double precision.

Compute distance with PBC.

void put_atoms_in_box ( int  ePBC,
const matrix  box,
gmx::ArrayRef< gmx::RVec x 
)

Put atoms inside the simulations box.

These routines puts ONE or ALL atoms in the box, not caring about charge groups! Also works for triclinic cells.

Parameters
[in]ePBCThe pbc type
[in]boxThe simulation box
[in,out]xThe coordinates of the atoms
void put_atoms_in_box_omp ( int  ePBC,
const matrix  box,
gmx::ArrayRef< gmx::RVec x,
int  nth 
)

Parallellizes put_atoms_in_box()

This wrapper function around put_atoms_in_box() with the ugly manual workload splitting is needed to avoid silently introducing multithreading in tools.

Parameters
[in]ePBCThe pbc type
[in]boxThe simulation box
[in,out]xThe coordinates of the atoms
[in]nthnumber of threads to be used in the given module
void put_atoms_in_compact_unitcell ( int  ePBC,
int  ecenter,
const matrix  box,
gmx::ArrayRef< gmx::RVec x 
)

Put atoms inside the unitcell.

This puts ALL atoms at the closest distance for the center of the box as calculated by calc_box_center. When ePBC=-1, the type of pbc is guessed from the box matrix.

Parameters
[in]ePBCThe pbc type
[in]ecenterThe pbc center type
[in]boxThe simulation box
[in,out]xThe coordinates of the atoms
void put_atoms_in_triclinic_unitcell ( int  ecenter,
const matrix  box,
gmx::ArrayRef< gmx::RVec x 
)

Put atoms inside triclinic box.

This puts ALL atoms in the triclinic unit cell, centered around the box center as calculated by calc_box_center.

Parameters
[in]ecenterThe pbc center type
[in]boxThe simulation box
[in,out]xThe coordinates of the atoms
void set_pbc ( t_pbc pbc,
int  ePBC,
const matrix  box 
)

Initiate the periodic boundary condition algorithms.

pbc_dx will not use pbc and return the normal difference vector when one or more of the diagonal elements of box are zero. When ePBC=-1, the type of pbc is guessed from the box matrix.

Parameters
[in,out]pbcThe pbc information structure
[in]ePBCThe PBC identifier
[in]boxThe box tensor
t_pbc* set_pbc_dd ( t_pbc pbc,
int  ePBC,
const ivec  domdecCells,
gmx_bool  bSingleDir,
const matrix  box 
)

Initiate the periodic boundary condition algorithms.

As set_pbc, but additionally sets that correct distances can be obtained using (combinations of) single box-vector shifts. Should be used with pbc_dx_aiuc. If domdecCells!=NULL pbc is not used for directions with dd->nc[i]==1 with bSingleDir==TRUE or with dd->nc[i]<=2 with bSingleDir==FALSE. Note that when no PBC is required only pbc->ePBC is set, the rest of the struct will be invalid.

Parameters
[in,out]pbcThe pbc information structure
[in]ePBCThe PBC identifier
[in]domdecCells3D integer vector describing the number of DD cells or nullptr if not using DD.
[in]bSingleDirTRUE if DD communicates only in one direction along dimensions
[in]boxThe box tensor
Returns
the pbc structure when pbc operations are required, NULL otherwise.