Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Macros | Functions
#include "gmxpre.h"
#include <cassert>
#include <cstdlib>
#include <cstring>
#include <algorithm>
#include <memory>
#include <string>
#include "gromacs/domdec/domdec.h"
#include "gromacs/domdec/domdec_network.h"
#include "gromacs/domdec/ga2la.h"
#include "gromacs/gmxlib/network.h"
#include "gromacs/math/vec.h"
#include "gromacs/mdlib/forcerec.h"
#include "gromacs/mdlib/gmx_omp_nthreads.h"
#include "gromacs/mdtypes/commrec.h"
#include "gromacs/mdtypes/inputrec.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/mdtypes/mdatom.h"
#include "gromacs/mdtypes/state.h"
#include "gromacs/pbcutil/mshift.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/topology/ifunc.h"
#include "gromacs/topology/mtop_util.h"
#include "gromacs/topology/topsort.h"
#include "gromacs/utility/cstringutil.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/logger.h"
#include "gromacs/utility/smalloc.h"
#include "gromacs/utility/strconvert.h"
#include "gromacs/utility/stringstream.h"
#include "gromacs/utility/stringutil.h"
#include "gromacs/utility/textwriter.h"
#include "domdec_constraints.h"
#include "domdec_internal.h"
#include "domdec_vsite.h"
#include "dump.h"
+ Include dependency graph for domdec_topology.cpp:

Description

This file defines functions used by the domdec module while managing the construction, use and error checking for topologies local to a DD rank.

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

Classes

struct  thread_work_t
 Struct for thread local work data for local topology generation. More...
 
struct  gmx_reverse_top_t
 Struct for the reverse topology: links bonded interactions to atomsx. More...
 

Macros

#define NITEM_DD_INIT_LOCAL_STATE   5
 The number of integer item in the local state, used for broadcasting of the state.
 

Functions

static int nral_rt (int ftype)
 Returns the number of atom entries for il in gmx_reverse_top_t.
 
static gmx_bool dd_check_ftype (int ftype, gmx_bool bBCheck, gmx_bool bConstr, gmx_bool bSettle)
 Return whether interactions of type ftype need to be assigned exactly once.
 
static std::string print_missing_interactions_mb (t_commrec *cr, const gmx_reverse_top_t *rt, const char *moltypename, const reverse_ilist_t *ril, int a_start, int a_end, int nat_mol, int nmol, const t_idef *idef)
 Help print error output when interactions are missing. More...
 
static void print_missing_interactions_atoms (const gmx::MDLogger &mdlog, t_commrec *cr, const gmx_mtop_t *mtop, const t_idef *idef)
 Help print error output when interactions are missing.
 
void dd_print_missing_interactions (const gmx::MDLogger &mdlog, t_commrec *cr, int local_count, const gmx_mtop_t *top_global, const gmx_localtop_t *top_local, const rvec *x, const matrix box)
 Print error output when interactions are missing.
 
static void global_atomnr_to_moltype_ind (const gmx_reverse_top_t *rt, int i_gl, int *mb, int *mt, int *mol, int *i_mol)
 Return global topology molecule information for global atom index i_gl.
 
static int getMaxNumExclusionsPerAtom (const t_blocka &excls)
 Returns the maximum number of exclusions per atom.
 
static int low_make_reverse_ilist (const InteractionLists &il_mt, const t_atom *atom, int *count, gmx_bool bConstr, gmx_bool bSettle, gmx_bool bBCheck, gmx::ArrayRef< const int > r_index, gmx::ArrayRef< int > r_il, gmx_bool bLinkToAllAtoms, gmx_bool bAssign)
 Run the reverse ilist generation and store it in r_il when bAssign = TRUE.
 
static int make_reverse_ilist (const InteractionLists &ilist, const t_atoms *atoms, gmx_bool bConstr, gmx_bool bSettle, gmx_bool bBCheck, gmx_bool bLinkToAllAtoms, reverse_ilist_t *ril_mt)
 Make the reverse ilist: a list of bonded interactions linked to atoms.
 
static gmx_reverse_top_t make_reverse_top (const gmx_mtop_t *mtop, gmx_bool bFE, gmx_bool bConstr, gmx_bool bSettle, gmx_bool bBCheck, int *nint)
 Generate the reverse topology.
 
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.
 
static void add_ifunc_for_vsites (t_iatom *tiatoms, const gmx_ga2la_t &ga2la, int nral, gmx_bool bHomeA, int a, int a_gl, int a_mol, const t_iatom *iatoms, t_ilist *il)
 Store a vsite interaction at the end of il. More...
 
static void add_ifunc (int nral, const t_iatom *tiatoms, t_ilist *il)
 Store a bonded interaction at the end of il.
 
static void add_posres (int mol, int a_mol, int numAtomsInMolecule, const gmx_molblock_t *molb, t_iatom *iatoms, const t_iparams *ip_in, t_idef *idef)
 Store a position restraint in idef and iatoms, complex because the parameters are different for each entry.
 
static void add_fbposres (int mol, int a_mol, int numAtomsInMolecule, const gmx_molblock_t *molb, t_iatom *iatoms, const t_iparams *ip_in, t_idef *idef)
 Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry.
 
static void add_vsite (const gmx_ga2la_t &ga2la, gmx::ArrayRef< const int > index, gmx::ArrayRef< const int > rtil, int ftype, int nral, gmx_bool bHomeA, int a, int a_gl, int a_mol, const t_iatom *iatoms, t_idef *idef)
 Store a virtual site interaction, complex because of PBC and recursion.
 
static real dd_dist2 (t_pbc *pbc_null, const rvec *x, const int i, int j)
 Returns the squared distance between atoms i and j.
 
static void combine_blocka (t_blocka *dest, gmx::ArrayRef< const thread_work_t > src)
 Append t_blocka block structures 1 to nsrc in src to *dest.
 
static void combine_idef (t_idef *dest, gmx::ArrayRef< const thread_work_t > src)
 Append t_idef structures 1 to nsrc in src to *dest.
 
static void check_assign_interactions_atom (int i, int i_gl, int mol, int i_mol, int numAtomsInMolecule, gmx::ArrayRef< const int > index, gmx::ArrayRef< const int > rtil, gmx_bool bInterMolInteractions, int ind_start, int ind_end, const gmx_domdec_t *dd, const gmx_domdec_zones_t *zones, const gmx_molblock_t *molb, gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B, real rc2, t_pbc *pbc_null, rvec *cg_cm, const t_iparams *ip_in, t_idef *idef, int iz, gmx_bool bBCheck, int *nbonded_local)
 Check and when available assign bonded interactions for local atom i.
 
static int make_bondeds_zone (gmx_domdec_t *dd, const gmx_domdec_zones_t *zones, const std::vector< gmx_molblock_t > &molb, gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B, real rc2, t_pbc *pbc_null, rvec *cg_cm, const t_iparams *ip_in, t_idef *idef, int izone, const gmx::Range< int > &atomRange)
 This function looks up and assigns bonded interactions for zone iz. More...
 
static void set_no_exclusions_zone (const gmx_domdec_zones_t *zones, int iz, t_blocka *lexcls)
 Set the exclusion data for i-zone iz for the case of no exclusions.
 
static void make_exclusions_zone (gmx_domdec_t *dd, gmx_domdec_zones_t *zones, const std::vector< gmx_moltype_t > &moltype, const int *cginfo, t_blocka *lexcls, int iz, int at_start, int at_end, const gmx::ArrayRef< const int > intermolecularExclusionGroup)
 Set the exclusion data for i-zone iz.
 
static void check_alloc_index (t_blocka *ba, int nindex_max)
 Ensure we have enough space in ba for nindex_max indices.
 
static void check_exclusions_alloc (const gmx_domdec_t *dd, const gmx_domdec_zones_t *zones, t_blocka *lexcls)
 Ensure that we have enough space for exclusion storate in lexcls.
 
static void finish_local_exclusions (gmx_domdec_t *dd, gmx_domdec_zones_t *zones, t_blocka *lexcls)
 Set the total count indexes for the local exclusions, needed by several functions.
 
static void clear_idef (t_idef *idef)
 Clear a t_idef struct.
 
static int make_local_bondeds_excls (gmx_domdec_t *dd, gmx_domdec_zones_t *zones, const gmx_mtop_t *mtop, const int *cginfo, gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B, real rc, t_pbc *pbc_null, rvec *cg_cm, t_idef *idef, t_blocka *lexcls, int *excl_count)
 Generate and store all required local bonded interactions in idef and local exclusions in lexcls.
 
void dd_make_local_top (gmx_domdec_t *dd, gmx_domdec_zones_t *zones, int npbcdim, matrix box, rvec cellsize_min, const ivec npulse, t_forcerec *fr, rvec *cgcm_or_x, const gmx_mtop_t &mtop, 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.
 
void dd_init_local_top (const gmx_mtop_t &top_global, gmx_localtop_t *top)
 Initialize local topology. More...
 
void dd_init_local_state (gmx_domdec_t *dd, const t_state *state_global, t_state *state_local)
 Construct local state.
 
static void check_link (t_blocka *link, int cg_gl, int cg_gl_j)
 Check if a link is stored in link between charge groups cg_gl and cg_gl_j and if not so, store a link.
 
t_blocka * makeBondedLinks (const gmx_mtop_t *mtop, cginfo_mb_t *cginfo_mb)
 Generate a list of links between atoms that are linked by bonded interactions. More...
 
static void update_max_bonded_distance (real r2, int ftype, int a1, int a2, bonded_distance_t *bd)
 Compare distance^2 r2 against the distance in bd and if larger store it along with ftype and atom indices a1 and a2.
 
static void bonded_cg_distance_mol (const gmx_moltype_t *molt, gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm, bonded_distance_t *bd_2b, bonded_distance_t *bd_mb)
 Set the distance, function type and atom indices for the longest distance between charge-groups of molecule type molt for two-body and multi-body bonded interactions.
 
static void bonded_distance_intermol (const InteractionLists &ilists_intermol, gmx_bool bBCheck, const rvec *x, int ePBC, const matrix box, bonded_distance_t *bd_2b, bonded_distance_t *bd_mb)
 Set the distance, function type and atom indices for the longest atom distance involved in intermolecular interactions for two-body and multi-body bonded interactions.
 
static bool moltypeHasVsite (const gmx_moltype_t &molt)
 Returns whether molt has at least one virtual site.
 
static void getWholeMoleculeCoordinates (const gmx_moltype_t *molt, const gmx_ffparams_t *ffparams, int ePBC, t_graph *graph, const matrix box, const rvec *x, rvec *xs)
 Returns coordinates not broken over PBC for a molecule.
 
void dd_bonded_cg_distance (const gmx::MDLogger &mdlog, 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.
 

Function Documentation

static void add_ifunc_for_vsites ( t_iatom *  tiatoms,
const gmx_ga2la_t ga2la,
int  nral,
gmx_bool  bHomeA,
int  a,
int  a_gl,
int  a_mol,
const t_iatom *  iatoms,
t_ilist *  il 
)
inlinestatic

Store a vsite interaction at the end of il.

This routine is very similar to add_ifunc, but vsites interactions have more work to do than other kinds of interactions, and the complex way nral (and thus vector contents) depends on ftype confuses static analysis tools unless we fuse the vsite atom-indexing organization code with the ifunc-adding code, so that they can see that nral is the same value.

void dd_init_local_top ( const gmx_mtop_t &  top_global,
gmx_localtop_t top 
)

Initialize local topology.

Parameters
[in]top_globalReference to global topology.
[in,out]topPointer to new local topology
static int make_bondeds_zone ( gmx_domdec_t *  dd,
const gmx_domdec_zones_t *  zones,
const std::vector< gmx_molblock_t > &  molb,
gmx_bool  bRCheckMB,
ivec  rcheck,
gmx_bool  bRCheck2B,
real  rc2,
t_pbc pbc_null,
rvec *  cg_cm,
const t_iparams *  ip_in,
t_idef *  idef,
int  izone,
const gmx::Range< int > &  atomRange 
)
static

This function looks up and assigns bonded interactions for zone iz.

With thread parallelizing each thread acts on a different atom range: at_start to at_end.

t_blocka* makeBondedLinks ( const gmx_mtop_t *  mtop,
cginfo_mb_t *  cginfo_mb 
)

Generate a list of links between atoms that are linked by bonded interactions.

Also stores whether atoms are linked in cginfo_mb.

static std::string print_missing_interactions_mb ( t_commrec *  cr,
const gmx_reverse_top_t rt,
const char *  moltypename,
const reverse_ilist_t *  ril,
int  a_start,
int  a_end,
int  nat_mol,
int  nmol,
const t_idef *  idef 
)
static

Help print error output when interactions are missing.

Note
This function needs to be called on all ranks (contains a global summation)