3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
45 #define GMX_LEFT 0 /* channel to the left processor */
46 #define GMX_RIGHT 1 /* channel to the right processor */
48 /* These are the good old ring communication routines */
50 extern void gmx_tx(const t_commrec
*cr
,int dir
,void *buf
,int bufsize
);
52 * Asynchronously sends bufsize bytes from the buffer pointed to by buf
53 * over the communication channel, identified by chan. The buffer becomes
54 * available after a successful call of gmx_tx_wait(dir).
57 extern void gmx_tx_wait(int dir
);
59 * Waits until the asynchronous send operation associated with chan has
60 * succeeded. This makes the buffer of the send operation available to
61 * the sending process.
64 extern void gmx_rx(const t_commrec
*cr
,int dir
,void *buf
,int bufsize
);
66 * Asynchronously receives bufsize bytes in the buffer pointed to by buf
67 * from communication channel identified by chan. The buffer becomes
68 * available after a successful call of gmx_rx_wait(chan).
71 extern void gmx_rx_wait(int dir
);
73 * Waits until the asynchronous receive operation, associated with chan,
74 * has succeeded. This makes the buffer of the receive operation
75 * available to the receiving process.
78 extern void gmx_left_right(int nnodes
,int nodeid
,
79 int *left
,int *right
);
80 /* Get left and right proc id. */
82 extern void gmx_tx_rx(const t_commrec
*cr
,
83 int send_dir
,void *send_buf
,int send_bufsize
,
84 int recv_dir
,void *recv_buf
,int recv_bufsize
);
85 /* Communicate simultaneously left and right */
87 extern void gmx_tx_rx_real(const t_commrec
*cr
,
88 int send_dir
,real
*send_buf
,int send_bufsize
,
89 int recv_dir
,real
*recv_buf
,int recv_bufsize
);
90 /* Communicate simultaneously left and right, reals only */
92 extern void gmx_wait(int dir_send
,int dir_recv
);
93 /* Wait for communication to finish */
95 extern void pd_move_f(const t_commrec
*cr
,rvec f
[],t_nrnb
*nrnb
);
96 /* Sum the forces over the nodes */
98 extern int *pd_cgindex(const t_commrec
*cr
);
100 extern int *pd_index(const t_commrec
*cr
);
102 extern int pd_shift(const t_commrec
*cr
);
104 extern int pd_bshift(const t_commrec
*cr
);
106 extern void pd_cg_range(const t_commrec
*cr
,int *cg0
,int *cg1
);
107 /* Get the range for the home charge groups */
109 extern void pd_at_range(const t_commrec
*cr
,int *at0
,int *at1
);
110 /* Get the range for the home particles */
112 extern gmx_localtop_t
*split_system(FILE *log
,
113 gmx_mtop_t
*mtop
,t_inputrec
*inputrec
,
115 /* Split the system over N processors. */
117 extern bool setup_parallel_vsites(t_idef
*idef
,t_commrec
*cr
,
118 t_comm_vsites
*vsitecomm
);
120 extern t_state
*partdec_init_local_state(t_commrec
*cr
,t_state
*state_global
);
121 /* Generate a local state struct from the global one */
124 pd_get_constraint_range(gmx_partdec_p_t pd
,int *start
,int *natoms
);
128 pd_constraints_nlocalatoms(gmx_partdec_p_t pd
);
130 /* Move x0 and also x1 if x1!=NULL */
132 pd_move_x_constraints(t_commrec
* cr
,