Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Functions
#include "gmxpre.h"
#include "domdec_constraints.h"
#include <cassert>
#include <algorithm>
#include <memory>
#include "gromacs/domdec/dlbtiming.h"
#include "gromacs/domdec/domdec.h"
#include "gromacs/domdec/domdec_struct.h"
#include "gromacs/domdec/ga2la.h"
#include "gromacs/domdec/hashedmap.h"
#include "gromacs/math/vec.h"
#include "gromacs/mdlib/constr.h"
#include "gromacs/mdlib/gmx_omp_nthreads.h"
#include "gromacs/mdtypes/commrec.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/pbcutil/ishift.h"
#include "gromacs/topology/ifunc.h"
#include "gromacs/topology/mtop_lookup.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/smalloc.h"
#include "domdec_internal.h"
#include "domdec_specatomcomm.h"
+ Include dependency graph for domdec_constraints.cpp:

Description

This file implements functions for domdec to use while managing inter-atomic constraints.

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

Classes

struct  gmx_domdec_constraints_t
 Struct used during constraint setup with domain decomposition. More...
 

Functions

void dd_move_x_constraints (gmx_domdec_t *dd, const 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.
 
gmx::ArrayRef< const int > dd_constraints_nlocalatoms (const gmx_domdec_t *dd)
 Returns the local atom count array for all constraints. More...
 
void dd_clear_local_constraint_indices (gmx_domdec_t *dd)
 Clears the local indices for the constraint communication setup.
 
static void walk_out (int con, int con_offset, int a, int offset, int nrec, gmx::ArrayRef< const int > ia1, gmx::ArrayRef< const int > ia2, const t_blocka *at2con, const gmx_ga2la_t &ga2la, gmx_bool bHomeConnect, gmx_domdec_constraints_t *dc, gmx_domdec_specat_comm_t *dcc, t_ilist *il_local, std::vector< int > *ireq)
 Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list.
 
static void atoms_to_settles (gmx_domdec_t *dd, const gmx_mtop_t *mtop, const int *cginfo, gmx::ArrayRef< const std::vector< int >> at2settle_mt, int cg_start, int cg_end, t_ilist *ils_local, std::vector< int > *ireq)
 Looks up SETTLE constraints for a range of charge-groups.
 
static void atoms_to_constraints (gmx_domdec_t *dd, const gmx_mtop_t *mtop, const int *cginfo, gmx::ArrayRef< const t_blocka > at2con_mt, int nrec, t_ilist *ilc_local, std::vector< int > *ireq)
 Looks up constraint for the local atoms.
 
int dd_make_local_constraints (gmx_domdec_t *dd, int at_start, const struct gmx_mtop_t *mtop, const int *cginfo, gmx::Constraints *constr, int nrec, t_ilist *il_local)
 Sets up communication and atom indices for all local+connected constraints.
 
void init_domdec_constraints (gmx_domdec_t *dd, const gmx_mtop_t *mtop)
 Initializes the data structures for constraint communication.
 

Function Documentation

gmx::ArrayRef<const int> dd_constraints_nlocalatoms ( const 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.

Note
When dd = nullptr, an empty reference is returned.