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
41 #include "types/commrec.h"
48 int * left_import_construct;
49 int left_import_nconstruct;
50 int * left_export_construct;
51 int left_export_nconstruct;
52 int * right_import_construct;
53 int right_import_nconstruct;
54 int * right_export_construct;
55 int right_export_nconstruct;
61 t_ilist ilist[F_NRE]; /* vsite ilists for this thread */
62 rvec fshift[SHIFTS]; /* fshift accumulation buffer */
63 matrix dxdf; /* virial dx*df accumulation buffer */
67 gmx_bool bHaveChargeGroups; /* Do we have charge groups? */
68 int n_intercg_vsite; /* The number of inter charge group vsites */
69 int nvsite_pbc_molt; /* The array size of vsite_pbc_molt */
70 int ***vsite_pbc_molt; /* The pbc atoms for intercg vsites */
71 int **vsite_pbc_loc; /* The local pbc atoms */
72 int *vsite_pbc_loc_nalloc; /* Sizes of vsite_pbc_loc */
73 gmx_bool bPDvsitecomm; /* Do we need vsite communication with PD? */
74 t_comm_vsites *vsitecomm; /* The PD vsite communication struct */
75 int nthreads; /* Number of threads used for vsites */
76 gmx_vsite_thread_t *tdata; /* Thread local vsites and work structs */
77 int *th_ind; /* Work array */
78 int th_ind_nalloc; /* Size of th_ind */
81 void construct_vsites(gmx_vsite_t *vsite,
84 t_iparams ip[], t_ilist ilist[],
85 int ePBC, gmx_bool bMolPBC, t_graph *graph,
86 t_commrec *cr, matrix box);
87 /* Create positions of vsite atoms based on surrounding atoms
88 * for the local system.
89 * If v is passed, the velocities of the vsites will be calculated
90 * as the new positions minus the old positions divided by dt,
91 * thus v should only be passed when the coordinates have been
92 * updated with a full time step.
93 * Note that velocitis of vsites are completely irrelevant
94 * for the integration, they are only useful for analysis.
97 void construct_vsites_mtop(gmx_vsite_t *vsite,
98 gmx_mtop_t *mtop, rvec x[]);
99 /* Create positions of vsite atoms based on surrounding atoms
100 * for the whole system.
101 * This function assumes that all molecules are whole.
104 void spread_vsite_f(gmx_vsite_t *vsite,
105 rvec x[], rvec f[], rvec *fshift,
106 gmx_bool VirCorr, matrix vir,
107 t_nrnb *nrnb, t_idef *idef,
108 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
110 /* Spread the force operating on the vsite atoms on the surrounding atoms.
111 * If fshift!=NULL also update the shift forces.
112 * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
113 * to vir. This correction is required when the virial is not calculated
114 * afterwards from the particle position and forces, but in a different way,
115 * as for instance for the PME mesh contribution.
118 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
119 gmx_bool bSerial_NoPBC);
120 /* Initialize the virtual site struct,
121 * returns NULL when there are no virtual sites.
122 * bSerial_NoPBC is to generate a simple vsite setup to be
123 * used only serial (no MPI or thread parallelization) and without pbc;
124 * this is useful for correction vsites of the initial configuration.
127 void split_vsites_over_threads(const t_ilist *ilist,
128 const t_mdatoms *mdatoms,
129 gmx_bool bLimitRange,
131 /* Divide the vsite work-load over the threads.
132 * Should be called at the end of the domain decomposition.
135 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
137 /* Set some vsite data for runs without domain decomposition.
138 * Should be called once after init_vsite, before calling other routines.