Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Macros | Functions
#include "gmxpre.h"
#include "partition.h"
#include "config.h"
#include <cassert>
#include <cstdio>
#include <algorithm>
#include "gromacs/domdec/collect.h"
#include "gromacs/domdec/dlb.h"
#include "gromacs/domdec/dlbtiming.h"
#include "gromacs/domdec/domdec_network.h"
#include "gromacs/domdec/ga2la.h"
#include "gromacs/domdec/localatomsetmanager.h"
#include "gromacs/domdec/localtopology.h"
#include "gromacs/domdec/localtopologychecker.h"
#include "gromacs/domdec/mdsetup.h"
#include "gromacs/domdec/nsgrid.h"
#include "gromacs/ewald/pme_pp.h"
#include "gromacs/gmxlib/network.h"
#include "gromacs/gmxlib/nrnb.h"
#include "gromacs/imd/imd.h"
#include "gromacs/math/functions.h"
#include "gromacs/math/vec.h"
#include "gromacs/mdlib/forcerec.h"
#include "gromacs/mdlib/gmx_omp_nthreads.h"
#include "gromacs/mdlib/mdatoms.h"
#include "gromacs/mdlib/vsite.h"
#include "gromacs/mdrunutility/mdmodulesnotifiers.h"
#include "gromacs/mdtypes/commrec.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/mdtypes/inputrec.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/mdtypes/mdatom.h"
#include "gromacs/mdtypes/nblist.h"
#include "gromacs/mdtypes/state.h"
#include "gromacs/nbnxm/nbnxm.h"
#include "gromacs/pulling/pull.h"
#include "gromacs/timing/wallcycle.h"
#include "gromacs/topology/mtop_util.h"
#include "gromacs/topology/topology.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/cstringutil.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/logger.h"
#include "gromacs/utility/real.h"
#include "gromacs/utility/strconvert.h"
#include "gromacs/utility/stringstream.h"
#include "gromacs/utility/stringutil.h"
#include "gromacs/utility/textwriter.h"
#include "box.h"
#include "cellsizes.h"
#include "distribute.h"
#include "domdec_constraints.h"
#include "domdec_internal.h"
#include "domdec_vsite.h"
#include "dump.h"
#include "redistribute.h"
#include "utility.h"
+ Include dependency graph for partition.cpp:

Description

This file defines functions for mdrun to call to make a new domain decomposition, and check it.

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

Classes

struct  dd_corners_t
 Domain corners for communication, a maximum of 4 i-zones see a j domain. More...
 

Macros

#define DD_PERF_LOSS_DLB_ON   0.02
 Turn on DLB when the load imbalance causes this amount of total loss. More...
 
#define DD_PERF_LOSS_WARN   0.05
 Warn about imbalance due to PP or PP/PME load imbalance at this loss.
 

Functions

static void print_ddzone (FILE *fp, int d, int i, int j, gmx_ddzone_t *zone)
 Debug helper printing a DD zone.
 
static void dd_move_cellx (gmx_domdec_t *dd, const gmx_ddbox_t *ddbox, rvec cell_ns_x0, rvec cell_ns_x1)
 Using the home grid size as input in cell_ns_x0 and cell_ns_x1 takes the extremes over all home and remote zones in the halo and returns the results in cell_ns_x0 and cell_ns_x1. Note: only used with the group cut-off scheme.
 
static void set_zones_numHomeAtoms (const gmx_domdec_t *dd)
 Sets the charge-group zones to be equal to the home zone.
 
static void restoreAtomGroups (gmx_domdec_t *dd, const t_state *state)
 Restore atom groups for the charge groups.
 
static void dd_set_atominfo (gmx::ArrayRef< const int > index_gl, int cg0, int cg1, t_forcerec *fr)
 Sets the atom info structures.
 
static void make_dd_indices (gmx_domdec_t *dd, const int atomStart)
 Makes the mappings between global and local atom indices during DD repartioning.
 
static void check_index_consistency (const gmx_domdec_t *dd, int natoms_sys, const char *where)
 Checks whether global and local atom indices are consistent.
 
static void clearDDStateIndices (gmx_domdec_t *dd, const bool keepLocalAtomIndices)
 Clear all DD global state indices.
 
static float dd_force_load (gmx_domdec_comm_t *comm)
 Return the duration of force calculations on this rank.
 
static void comm_dd_ns_cell_sizes (gmx_domdec_t *dd, gmx_ddbox_t *ddbox, rvec cell_ns_x0, rvec cell_ns_x1, int64_t step)
 Runs cell size checks and communicates the boundaries.
 
static void get_load_distribution (gmx_domdec_t *dd, gmx_wallcycle *wcycle)
 Compute and communicate to determine the load distribution across PP ranks.
 
static float dd_force_load_fraction (gmx_domdec_t *dd)
 Return the relative performance loss on the total run time due to the force calculation load imbalance.
 
static float dd_force_imb_perf_loss (gmx_domdec_t *dd)
 Return the relative performance loss on the total run time due to the force calculation load imbalance.
 
static void print_dd_load_av (FILE *fplog, gmx_domdec_t *dd)
 Print load-balance report e.g. at the end of a run.
 
static float dd_vol_min (gmx_domdec_t *dd)
 Return the minimum communication volume.
 
static int dd_load_flags (gmx_domdec_t *dd)
 Return the DD load flags.
 
static float dd_f_imbal (gmx_domdec_t *dd)
 Return the reported load imbalance in force calculations.
 
static std::string dd_print_load (gmx_domdec_t *dd, int64_t step)
 Returns DD load balance report.
 
static void dd_print_load_verbose (gmx_domdec_t *dd)
 Prints DD load balance report in mdrun verbose mode.
 
static void turn_on_dlb (const gmx::MDLogger &mdlog, gmx_domdec_t *dd, int64_t step)
 Turns on dynamic load balancing if possible and needed.
 
static void turn_off_dlb (const gmx::MDLogger &mdlog, gmx_domdec_t *dd, int64_t step)
 Turns off dynamic load balancing (but leave it able to turn back on).
 
static void turn_off_dlb_forever (const gmx::MDLogger &mdlog, gmx_domdec_t *dd, int64_t step)
 Turns off dynamic load balancing permanently.
 
void set_dd_dlb_max_cutoff (t_commrec *cr, real cutoff)
 Limit DLB to preserve the option of returning to the current cut-off. More...
 
static void merge_cg_buffers (int ncell, gmx_domdec_comm_dim_t *cd, int pulse, int *ncg_cell, gmx::ArrayRef< int > index_gl, const int *recv_i, gmx::ArrayRef< gmx::RVec > x, gmx::ArrayRef< const gmx::RVec > recv_vr, gmx::ArrayRef< gmx::AtomInfoWithinMoleculeBlock > atomInfoForEachMoleculeBlock, gmx::ArrayRef< int64_t > atomInfo)
 Merge atom buffers.
 
static void make_cell2at_index (gmx_domdec_comm_dim_t *cd, int nzone, int atomGroupStart)
 Makes a range partitioning for the atom groups wthin a cell.
 
static bool missing_link (const gmx::ListOfLists< int > &link, const int globalAtomIndex, const gmx_ga2la_t &ga2la)
 Returns whether a link is missing.
 
static void set_dd_corners (const gmx_domdec_t *dd, int dim0, int dim1, int dim2, gmx_bool bDistMB, dd_corners_t *c)
 Determine the corners of the domain(s) we are communicating with.
 
static void get_zone_pulse_groups (gmx_domdec_t *dd, int zonei, int zone, int cg0, int cg1, gmx::ArrayRef< const int > globalAtomGroupIndices, int dim, int dim_ind, int dim0, int dim1, int dim2, real r_comm2, real r_bcomm2, matrix box, bool distanceIsTriclinic, rvec *normal, real skew_fac2_d, real skew_fac_01, rvec *v_d, rvec *v_0, rvec *v_1, const dd_corners_t *c, const rvec sf2_round, gmx_bool bDistBonded, gmx_bool bBondComm, gmx_bool bDist2B, gmx_bool bDistMB, gmx::ArrayRef< const gmx::RVec > coordinates, gmx::ArrayRef< const int64_t > atomInfo, std::vector< int > *localAtomGroups, dd_comm_setup_work_t *work)
 Add the atom groups and coordinates we need to send in this pulse from this zone to localAtomGroups and work.
 
static void clearCommSetupData (dd_comm_setup_work_t *work)
 Clear data.
 
static void setup_dd_communication (gmx_domdec_t *dd, matrix box, gmx_ddbox_t *ddbox, t_forcerec *fr, t_state *state)
 Prepare DD communication.
 
static void set_cg_boundaries (gmx_domdec_zones_t *zones)
 Set boundaries for the charge group range.
 
static void set_zones_size (gmx_domdec_t *dd, matrix box, const gmx_ddbox_t *ddbox, int zone_start, int zone_end, int numMovedChargeGroupsInHomeZone)
 Set zone dimensions for zones zone_start to zone_end-1. More...
 
template<typename T >
static void orderVector (gmx::ArrayRef< const gmx_cgsort_t > sort, gmx::ArrayRef< T > dataToSort, gmx::ArrayRef< T > sortBuffer)
 Order data in dataToSort according to sort. More...
 
template<typename T >
static void orderVector (gmx::ArrayRef< const gmx_cgsort_t > sort, gmx::ArrayRef< T > vectorToSort, std::vector< T > *workVector)
 Order data in dataToSort according to sort. More...
 
static void dd_sort_order_nbnxn (const t_forcerec *fr, std::vector< gmx_cgsort_t > *sort)
 Returns the sorting order for atoms based on the nbnxn grid order in sort.
 
static void dd_sort_state (gmx_domdec_t *dd, t_forcerec *fr, t_state *state)
 Returns the sorting state for DD.
 
static void add_dd_statistics (gmx_domdec_t *dd)
 Accumulates load statistics.
 
void reset_dd_statistics_counters (gmx_domdec_t *dd)
 Reset all the statistics and counters for total run counting.
 
bool gmx::check_grid_jump (int64_t step, const gmx_domdec_t *dd, real cutoff, const gmx_ddbox_t *ddbox, bool bFatal)
 Check whether the DD grid has moved too far for correctness.
 
void gmx::print_dd_statistics (const t_commrec *cr, const t_inputrec &inputrec, FILE *fplog)
 Print statistics for domain decomposition communication.
 
void gmx::dd_partition_system (FILE *fplog, const gmx::MDLogger &mdlog, int64_t step, const t_commrec *cr, bool bMainState, t_state *state_global, const gmx_mtop_t &top_global, const t_inputrec &inputrec, const MDModulesNotifiers &mdModulesNotifiers, gmx::ImdSession *imdSession, pull_t *pull_work, t_state *state_local, gmx::ForceBuffers *f, gmx::MDAtoms *mdAtoms, gmx_localtop_t *top_local, t_forcerec *fr, gmx::VirtualSitesHandler *vsite, gmx::Constraints *constr, t_nrnb *nrnb, gmx_wallcycle *wcycle, bool bVerbose)
 TODO Remove fplog when group scheme and charge groups are gone. More...
 

Macro Definition Documentation

#define DD_PERF_LOSS_DLB_ON   0.02

Turn on DLB when the load imbalance causes this amount of total loss.

There is a bit of overhead with DLB and it's difficult to achieve a load imbalance of less than 2% with DLB.

Function Documentation

template<typename T >
static void orderVector ( gmx::ArrayRef< const gmx_cgsort_t >  sort,
gmx::ArrayRef< T >  dataToSort,
gmx::ArrayRef< T >  sortBuffer 
)
static

Order data in dataToSort according to sort.

Note: both buffers should have at least sort.size() elements.

template<typename T >
static void orderVector ( gmx::ArrayRef< const gmx_cgsort_t >  sort,
gmx::ArrayRef< T >  vectorToSort,
std::vector< T > *  workVector 
)
static

Order data in dataToSort according to sort.

Note: vectorToSort should have at least sort.size() elements, workVector is resized when it is too small.

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.

static void set_zones_size ( gmx_domdec_t *  dd,
matrix  box,
const gmx_ddbox_t *  ddbox,
int  zone_start,
int  zone_end,
int  numMovedChargeGroupsInHomeZone 
)
static

Set zone dimensions for zones zone_start to zone_end-1.

Also sets the atom density for the home zone when zone_start=0. For this numMovedChargeGroupsInHomeZone needs to be passed to tell how many charge groups will move but are still part of the current range.

Todo:
When converting domdec to use proper classes, all these variables should be private and a method should return the correct count depending on an internal state.
Parameters
[in,out]ddThe domain decomposition struct
[in]boxThe box
[in]ddboxThe domain decomposition box struct
[in]zone_startThe start of the zone range to set sizes for
[in]zone_endThe end of the zone range to set sizes for
[in]numMovedChargeGroupsInHomeZoneThe number of charge groups in the home zone that should moved but are still present in dd->comm->zones.cg_range