Gromacs
2018.8
|
#include <stdio.h>
#include <vector>
#include "gromacs/gmxlib/nrnb.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/mdlib/vsite.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/mdtypes/mdatom.h"
#include "gromacs/timing/wallcycle.h"
#include "gromacs/topology/block.h"
#include "gromacs/topology/idef.h"
#include "gromacs/topology/topology.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/basedefinitions.h"
#include "gromacs/utility/real.h"
This file declares functions for mdrun to call to manage the details of its domain decomposition.
Classes | |
struct | DomdecOptions |
Structure containing all (command line) options for the domain decomposition. More... | |
Enumerations | |
enum | DdRankOrder { DdRankOrder::select, DdRankOrder::interleave, DdRankOrder::pp_pme, DdRankOrder::cartesian, DdRankOrder::nr } |
The options for the domain decomposition MPI task ordering. More... | |
enum | DlbOption { DlbOption::select, DlbOption::turnOnWhenUseful, DlbOption::no, DlbOption::yes, DlbOption::nr } |
The options for the dynamic load balancing. More... | |
enum | { ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr } |
Cycle counter indices used internally in the domain decomposition. | |
Functions | |
int | ddglatnr (const gmx_domdec_t *dd, int i) |
Returns the global topology atom number belonging to local atom index i. More... | |
t_block * | dd_charge_groups_global (struct gmx_domdec_t *dd) |
Return a block struct for the charge groups of the whole system. | |
void | dd_store_state (struct gmx_domdec_t *dd, t_state *state) |
Store the global cg indices of the home cgs in state,. More... | |
struct gmx_domdec_zones_t * | domdec_zones (struct gmx_domdec_t *dd) |
Returns a pointer to the gmx_domdec_zones_t struct. | |
void | dd_get_ns_ranges (const gmx_domdec_t *dd, int icg, int *jcg0, int *jcg1, ivec shift0, ivec shift1) |
Sets the j-charge-group range for i-charge-group icg . | |
int | dd_natoms_mdatoms (const gmx_domdec_t *dd) |
Returns the atom range in the local state for atoms that need to be present in mdatoms. | |
int | dd_natoms_vsite (const gmx_domdec_t *dd) |
Returns the atom range in the local state for atoms involved in virtual sites. | |
void | dd_get_constraint_range (const gmx_domdec_t *dd, int *at_start, int *at_end) |
Sets the atom range for atom in the local state for atoms received in constraints communication. | |
void | get_pme_nnodes (const struct gmx_domdec_t *dd, int *npmenodes_x, int *npmenodes_y) |
Get the number of PME nodes along x and y, can be called with dd=NULL. | |
std::vector< int > | get_pme_ddranks (t_commrec *cr, int pmenodeid) |
Returns the set of DD ranks that communicate with pme node cr->nodeid. | |
int | dd_pme_maxshift_x (const gmx_domdec_t *dd) |
Returns the maximum shift for coordinate communication in PME, dim x. | |
int | dd_pme_maxshift_y (const gmx_domdec_t *dd) |
Returns the maximum shift for coordinate communication in PME, dim y. | |
gmx_domdec_t * | init_domain_decomposition (FILE *fplog, t_commrec *cr, const DomdecOptions &options, const MdrunOptions &mdrunOptions, const gmx_mtop_t *mtop, const t_inputrec *ir, const matrix box, const rvec *xGlobal, gmx_ddbox_t *ddbox, int *npme_x, int *npme_y) |
Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct. | |
void | dd_init_bondeds (FILE *fplog, gmx_domdec_t *dd, const gmx_mtop_t *mtop, const gmx_vsite_t *vsite, const t_inputrec *ir, gmx_bool bBCheck, cginfo_mb_t *cginfo_mb) |
Initialize data structures for bonded interactions. | |
gmx_bool | dd_bonded_molpbc (const gmx_domdec_t *dd, int ePBC) |
Returns if we need to do pbc for calculating bonded interactions. | |
gmx_bool | change_dd_cutoff (struct t_commrec *cr, t_state *state, const t_inputrec *ir, real cutoff_req) |
Change the DD non-bonded communication cut-off. More... | |
void | set_dd_dlb_max_cutoff (struct t_commrec *cr, real cutoff) |
Limit DLB to preserve the option of returning to the current cut-off. More... | |
gmx_bool | dd_dlb_is_on (const struct gmx_domdec_t *dd) |
Return if we are currently using dynamic load balancing. | |
gmx_bool | dd_dlb_is_locked (const struct gmx_domdec_t *dd) |
Return if the DLB lock is set. | |
void | dd_dlb_lock (struct gmx_domdec_t *dd) |
Set a lock such that with DLB=auto DLB cannot get turned on. | |
void | dd_dlb_unlock (struct gmx_domdec_t *dd) |
Clear a lock such that with DLB=auto DLB may get turned on later. | |
void | dd_setup_dlb_resource_sharing (t_commrec *cr, int gpu_id) |
Set up communication for averaging GPU wait times over domains. More... | |
void | dd_collect_vec (struct gmx_domdec_t *dd, const t_state *state_local, gmx::ArrayRef< const gmx::RVec > lv, gmx::ArrayRef< gmx::RVec > v) |
Collects local rvec arrays lv to v on the master rank. | |
void | dd_collect_state (struct gmx_domdec_t *dd, const t_state *state_local, t_state *state) |
Collects the local state state_local to state on the master rank. | |
void | dd_cycles_add (const gmx_domdec_t *dd, float cycles, int ddCycl) |
Add the wallcycle count to the DD counter. | |
void | dd_force_flop_start (struct gmx_domdec_t *dd, t_nrnb *nrnb) |
Start the force flop count. | |
void | dd_force_flop_stop (struct gmx_domdec_t *dd, t_nrnb *nrnb) |
Stop the force flop count. | |
float | dd_pme_f_ratio (struct gmx_domdec_t *dd) |
Return the PME/PP force load ratio, or -1 if nothing was measured. More... | |
void | dd_move_x (struct gmx_domdec_t *dd, matrix box, rvec x[]) |
Communicate the coordinates to the neighboring cells and do pbc. | |
void | dd_move_f (struct gmx_domdec_t *dd, rvec f[], rvec *fshift) |
Sum the forces over the neighboring cells. More... | |
void | dd_atom_spread_real (struct gmx_domdec_t *dd, real v[]) |
Communicate a real for each atom to the neighboring cells. | |
void | dd_atom_sum_real (struct gmx_domdec_t *dd, real v[]) |
Sum the contributions to a real for each atom over the neighboring cells. | |
void | dd_partition_system (FILE *fplog, gmx_int64_t step, t_commrec *cr, gmx_bool bMasterState, int nstglobalcomm, t_state *state_global, const gmx_mtop_t *top_global, const t_inputrec *ir, t_state *state_local, PaddedRVecVector *f, gmx::MDAtoms *mdatoms, gmx_localtop_t *top_local, t_forcerec *fr, gmx_vsite_t *vsite, struct gmx_constr *constr, t_nrnb *nrnb, gmx_wallcycle_t wcycle, gmx_bool bVerbose) |
Partition the system over the nodes. More... | |
void | reset_dd_statistics_counters (struct gmx_domdec_t *dd) |
Reset all the statistics and counters for total run counting. | |
void | print_dd_statistics (struct t_commrec *cr, const t_inputrec *ir, FILE *fplog) |
Print statistics for domain decomposition communication. | |
void | dd_move_f_vsites (struct gmx_domdec_t *dd, rvec *f, rvec *fshift) |
Communicates the virtual site forces, reduces the shift forces when fshift != NULL. | |
void | dd_clear_f_vsites (struct gmx_domdec_t *dd, rvec *f) |
Clears the forces for virtual sites. | |
void | dd_move_x_constraints (struct gmx_domdec_t *dd, matrix box, rvec *x0, rvec *x1, gmx_bool bX1IsCoord) |
Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1. | |
void | dd_move_x_vsites (struct gmx_domdec_t *dd, const matrix box, rvec *x) |
Communicates the coordinates involved in virtual sites. | |
int * | dd_constraints_nlocalatoms (struct gmx_domdec_t *dd) |
Returns the local atom count array for all constraints. More... | |
void | dd_print_missing_interactions (FILE *fplog, struct t_commrec *cr, int local_count, const gmx_mtop_t *top_global, const gmx_localtop_t *top_local, t_state *state_local) |
Print error output when interactions are missing. | |
void | dd_make_reverse_top (FILE *fplog, gmx_domdec_t *dd, const gmx_mtop_t *mtop, const gmx_vsite_t *vsite, const t_inputrec *ir, gmx_bool bBCheck) |
Generate and store the reverse topology. | |
void | dd_make_local_cgs (struct gmx_domdec_t *dd, t_block *lcgs) |
Store the local charge group index in lcgs . | |
void | dd_make_local_top (struct gmx_domdec_t *dd, struct gmx_domdec_zones_t *zones, int npbcdim, matrix box, rvec cellsize_min, ivec npulse, t_forcerec *fr, rvec *cgcm_or_x, gmx_vsite_t *vsite, const gmx_mtop_t *top, gmx_localtop_t *ltop) |
Generate the local topology and virtual site data. | |
void | dd_sort_local_top (gmx_domdec_t *dd, const t_mdatoms *mdatoms, gmx_localtop_t *ltop) |
Sort ltop->ilist when we are doing free energy. | |
gmx_localtop_t * | dd_init_local_top (const gmx_mtop_t *top_global) |
Construct local topology. | |
void | dd_init_local_state (struct gmx_domdec_t *dd, t_state *state_global, t_state *local_state) |
Construct local state. | |
t_blocka * | make_charge_group_links (const gmx_mtop_t *mtop, gmx_domdec_t *dd, cginfo_mb_t *cginfo_mb) |
Generate a list of links between charge groups that are linked by bonded interactions. | |
void | dd_bonded_cg_distance (FILE *fplog, const gmx_mtop_t *mtop, const t_inputrec *ir, const rvec *x, const matrix box, gmx_bool bBCheck, real *r_2b, real *r_mb) |
Calculate the maximum distance involved in 2-body and multi-body bonded interactions. | |
void | write_dd_pdb (const char *fn, gmx_int64_t step, const char *title, const gmx_mtop_t *mtop, t_commrec *cr, int natoms, rvec x[], matrix box) |
Dump a pdb file with the current DD home + communicated atoms. More... | |
real | comm_box_frac (const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox) |
Returns the volume fraction of the system that is communicated. | |
real | dd_choose_grid (FILE *fplog, t_commrec *cr, gmx_domdec_t *dd, const t_inputrec *ir, const gmx_mtop_t *mtop, const matrix box, const gmx_ddbox_t *ddbox, int nPmeRanks, gmx_bool bDynLoadBal, real dlb_scale, real cellsize_limit, real cutoff_dd, gmx_bool bInterCGBondeds) |
Determines the optimal DD cell setup dd->nc and possibly npmenodes for the system. More... | |
void | set_ddbox (gmx_domdec_t *dd, gmx_bool bMasterState, t_commrec *cr_sum, const t_inputrec *ir, const matrix box, gmx_bool bCalcUnboundedSize, const t_block *cgs, const rvec *x, gmx_ddbox_t *ddbox) |
Set the box and PBC data in ddbox . | |
void | set_ddbox_cr (t_commrec *cr, const ivec *dd_nc, const t_inputrec *ir, const matrix box, const t_block *cgs, const rvec *x, gmx_ddbox_t *ddbox) |
Set the box and PBC data in ddbox . | |
|
strong |
The options for the domain decomposition MPI task ordering.
|
strong |
gmx_bool change_dd_cutoff | ( | struct t_commrec * | cr, |
t_state * | state, | ||
const t_inputrec * | ir, | ||
real | cutoff_req | ||
) |
Change the DD non-bonded communication cut-off.
This could fail when trying to increase the cut-off, then FALSE will be returned and the cut-off is not modified.
real dd_choose_grid | ( | FILE * | fplog, |
t_commrec * | cr, | ||
gmx_domdec_t * | dd, | ||
const t_inputrec * | ir, | ||
const gmx_mtop_t * | mtop, | ||
const matrix | box, | ||
const gmx_ddbox_t * | ddbox, | ||
int | nPmeRanks, | ||
gmx_bool | bDynLoadBal, | ||
real | dlb_scale, | ||
real | cellsize_limit, | ||
real | cutoff_dd, | ||
gmx_bool | bInterCGBondeds | ||
) |
Determines the optimal DD cell setup dd->nc and possibly npmenodes for the system.
On the master node returns the actual cellsize limit used.
int* dd_constraints_nlocalatoms | ( | struct gmx_domdec_t * | dd | ) |
Returns the local atom count array for all constraints.
The local atom count for a constraint, possible values 2/1/0, is needed to avoid not/double-counting contributions linked to the Lagrange multiplier, such as the virial and free-energy derivatives.
void dd_move_f | ( | struct gmx_domdec_t * | dd, |
rvec | f[], | ||
rvec * | fshift | ||
) |
Sum the forces over the neighboring cells.
When fshift!=NULL the shift forces are updated to obtain the correct virial from the single sum including f.
void dd_partition_system | ( | FILE * | fplog, |
gmx_int64_t | step, | ||
t_commrec * | cr, | ||
gmx_bool | bMasterState, | ||
int | nstglobalcomm, | ||
t_state * | state_global, | ||
const gmx_mtop_t * | top_global, | ||
const t_inputrec * | ir, | ||
t_state * | state_local, | ||
PaddedRVecVector * | f, | ||
gmx::MDAtoms * | mdatoms, | ||
gmx_localtop_t * | top_local, | ||
t_forcerec * | fr, | ||
gmx_vsite_t * | vsite, | ||
struct gmx_constr * | constr, | ||
t_nrnb * | nrnb, | ||
gmx_wallcycle_t | wcycle, | ||
gmx_bool | bVerbose | ||
) |
Partition the system over the nodes.
step is only used for printing error messages. If bMasterState==TRUE then state_global from the master node is used, else state_local is redistributed between the nodes. When f!=NULL, *f will be reallocated to the size of state_local.
float dd_pme_f_ratio | ( | struct gmx_domdec_t * | dd | ) |
Return the PME/PP force load ratio, or -1 if nothing was measured.
Should only be called on the DD master node.
void dd_setup_dlb_resource_sharing | ( | t_commrec * | cr, |
int | gpu_id | ||
) |
Set up communication for averaging GPU wait times over domains.
When domains (PP MPI ranks) share a GPU, the individual GPU wait times are meaningless, as it depends on the order in which tasks on the same GPU finish. Therefore there wait times need to be averaged over the ranks sharing the same GPU. This function sets up the communication for that.
void dd_store_state | ( | struct gmx_domdec_t * | dd, |
t_state * | state | ||
) |
Store the global cg indices of the home cgs in state,.
This means it can be reset, even after a new DD partitioning.
int ddglatnr | ( | const gmx_domdec_t * | dd, |
int | i | ||
) |
Returns the global topology atom number belonging to local atom index i.
This function is intended for writing ASCII output and returns atom numbers starting at 1. When dd=NULL returns i+1.
void set_dd_dlb_max_cutoff | ( | struct t_commrec * | cr, |
real | cutoff | ||
) |
Limit DLB to preserve the option of returning to the current cut-off.
Domain boundary changes due to the DD dynamic load balancing can limit the cut-off distance that can be set in change_dd_cutoff. This function sets/changes the DLB limit such that using the passed (pair-list) cut-off should still be possible after subsequently setting a shorter cut-off with change_dd_cutoff.
void write_dd_pdb | ( | const char * | fn, |
gmx_int64_t | step, | ||
const char * | title, | ||
const gmx_mtop_t * | mtop, | ||
t_commrec * | cr, | ||
int | natoms, | ||
rvec | x[], | ||
matrix | box | ||
) |
Dump a pdb file with the current DD home + communicated atoms.
When natoms=-1, dump all known atoms.