|
Gromacs
2026.0-dev-20251119-5f0a571d
|
#include "gmxpre.h"#include "domdec_constraints.h"#include <cassert>#include <cstdio>#include <algorithm>#include <array>#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/mdlib/constr.h"#include "gromacs/mdlib/gmx_omp_nthreads.h"#include "gromacs/mdtypes/atominfo.h"#include "gromacs/mdtypes/commrec.h"#include "gromacs/pbcutil/ishift.h"#include "gromacs/topology/idef.h"#include "gromacs/topology/mtop_lookup.h"#include "gromacs/topology/topology.h"#include "gromacs/utility/basedefinitions.h"#include "gromacs/utility/exceptions.h"#include "gromacs/utility/fatalerror.h"#include "gromacs/utility/gmxassert.h"#include "gromacs/utility/listoflists.h"#include "gromacs/utility/vec.h"#include "gromacs/utility/vectypes.h"#include "domdec_internal.h"#include "domdec_specatomcomm.h"
Include dependency graph for domdec_constraints.cpp:This file implements functions for domdec to use while managing inter-atomic constraints.
Functions | |
| void | dd_move_x_constraints (gmx_domdec_t *dd, const matrix box, gmx::ArrayRef< gmx::RVec > x0, gmx::ArrayRef< gmx::RVec > x1, 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 ListOfLists< int > &at2con, const gmx_ga2la_t &ga2la, gmx_bool bHomeConnect, gmx_domdec_constraints_t *dc, gmx_domdec_specat_comm_t *dcc, InteractionList *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, gmx::ArrayRef< const int32_t > atomInfo, gmx::ArrayRef< const std::vector< int >> at2settle_mt, int cg_start, int cg_end, InteractionList *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, gmx::ArrayRef< const int > atomInfo, gmx::ArrayRef< const ListOfLists< int >> at2con_mt, int nrec, InteractionList *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, gmx::ArrayRef< const int32_t > atomInfo, gmx::Constraints *constr, int nrec, gmx::EnumerationArray< InteractionFunction, InteractionList > &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. | |
| 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.
dd = nullptr, an empty reference is returned.
1.8.5