#include "gmxpre.h"
#include "domdec_specatomcomm.h"
#include <assert.h>
#include <algorithm>
#include "gromacs/domdec/domdec.h"
#include "gromacs/domdec/domdec_network.h"
#include "gromacs/legacyheaders/gmx_ga2la.h"
#include "gromacs/legacyheaders/gmx_hash.h"
#include "gromacs/legacyheaders/gmx_omp_nthreads.h"
#include "gromacs/legacyheaders/macros.h"
#include "gromacs/legacyheaders/types/commrec.h"
#include "gromacs/math/vec.h"
#include "gromacs/pbcutil/ishift.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/smalloc.h"
This file implements functions for domdec to use while managing inter-atomic constraints.
- Author
- Berk Hess hess@.nosp@m.kth..nosp@m.se
|
void | dd_move_f_specat (gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac, rvec *f, rvec *fshift) |
| Communicates the force for special atoms, the shift forces are reduced with fshift != NULL.
|
|
void | dd_move_x_specat (gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac, matrix box, rvec *x0, rvec *x1, gmx_bool bX1IsCoord) |
| Communicates the coordinates for special atoms. More...
|
|
int | setup_specat_communication (gmx_domdec_t *dd, ind_req_t *ireq, gmx_domdec_specat_comm_t *spac, gmx_hash_t ga2la_specat, int at_start, int vbuf_fac, const char *specat_type, const char *add_err) |
| Sets up the communication for special atoms. More...
|
|
gmx_domdec_specat_comm_t * | specat_comm_init (int nthread) |
| Initialize a special communication struct.
|
|
Communicates the coordinates for special atoms.
- Parameters
-
[in] | dd | Domain decomposition struct |
[in] | spac | Special atom communication struct |
[in] | box | Box, used for pbc |
[in,out] | x0 | Vector to communicate |
[in,out] | x1 | Vector to communicate, when != NULL |
[in] | bX1IsCoord | Tells is x1 is a coordinate vector (needs pbc) |
int setup_specat_communication |
( |
gmx_domdec_t * |
dd, |
|
|
ind_req_t * |
ireq, |
|
|
gmx_domdec_specat_comm_t * |
spac, |
|
|
gmx_hash_t |
ga2la_specat, |
|
|
int |
at_start, |
|
|
int |
vbuf_fac, |
|
|
const char * |
specat_type, |
|
|
const char * |
add_err |
|
) |
| |
Sets up the communication for special atoms.
- Parameters
-
[in] | dd | Domain decomposition struct |
[in] | ireq | List of requested atom indices |
[in,out] | spac | Special atom communication struct |
[out] | ga2la_specat | Global to local special atom index |
[in] | at_start | Index in local state where to start storing communicated atoms |
[in] | vbuf_fac | Buffer factor, 1 or 2 for communicating 1 or 2 vectors |
[in] | specat_type | Name of the special atom, used for error message |
[in] | add_err | Text to add at the end of error message when atoms can't be found |