Run regressiontests from build
[gromacs.git] / include / vsite.h
blob50edc3d405caccd346f4a5512c444c2be2b88944
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #ifndef _vsite_h
40 #define _vsite_h
41 #include <stdio.h>
42 #include "visibility.h"
43 #include "typedefs.h"
44 #include "types/commrec.h"
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
50 typedef struct {
51 int * left_import_construct;
52 int left_import_nconstruct;
53 int * left_export_construct;
54 int left_export_nconstruct;
55 int * right_import_construct;
56 int right_import_nconstruct;
57 int * right_export_construct;
58 int right_export_nconstruct;
59 rvec * send_buf;
60 rvec * recv_buf;
61 } t_comm_vsites;
63 typedef struct {
64 t_ilist ilist[F_NRE]; /* vsite ilists for this thread */
65 rvec fshift[SHIFTS]; /* fshift accumulation buffer */
66 matrix dxdf; /* virial dx*df accumulation buffer */
67 } gmx_vsite_thread_t;
69 typedef struct {
70 gmx_bool bHaveChargeGroups; /* Do we have charge groups? */
71 int n_intercg_vsite; /* The number of inter charge group vsites */
72 int nvsite_pbc_molt; /* The array size of vsite_pbc_molt */
73 int ***vsite_pbc_molt; /* The pbc atoms for intercg vsites */
74 int **vsite_pbc_loc; /* The local pbc atoms */
75 int *vsite_pbc_loc_nalloc; /* Sizes of vsite_pbc_loc */
76 gmx_bool bPDvsitecomm; /* Do we need vsite communication with PD? */
77 t_comm_vsites *vsitecomm; /* The PD vsite communication struct */
78 int nthreads; /* Number of threads used for vsites */
79 gmx_vsite_thread_t *tdata; /* Thread local vsites and work structs */
80 int *th_ind; /* Work array */
81 int th_ind_nalloc; /* Size of th_ind */
82 } gmx_vsite_t;
84 GMX_LIBMD_EXPORT
85 void construct_vsites(FILE *log,gmx_vsite_t *vsite,
86 rvec x[],t_nrnb *nrnb,
87 real dt,rvec v[],
88 t_iparams ip[],t_ilist ilist[],
89 int ePBC,gmx_bool bMolPBC,t_graph *graph,
90 t_commrec *cr,matrix box);
91 /* Create positions of vsite atoms based on surrounding atoms
92 * for the local system.
93 * If v is passed, the velocities of the vsites will be calculated
94 * as the new positions minus the old positions divided by dt,
95 * thus v should only be passed when the coordinates have been
96 * updated with a full time step.
97 * Note that velocitis of vsites are completely irrelevant
98 * for the integration, they are only useful for analysis.
101 GMX_LIBMD_EXPORT
102 void construct_vsites_mtop(FILE *log,gmx_vsite_t *vsite,
103 gmx_mtop_t *mtop,rvec x[]);
104 /* Create positions of vsite atoms based on surrounding atoms
105 * for the whole system.
106 * This function assumes that all molecules are whole.
109 void spread_vsite_f(FILE *log,gmx_vsite_t *vsite,
110 rvec x[],rvec f[],rvec *fshift,
111 gmx_bool VirCorr,matrix vir,
112 t_nrnb *nrnb,t_idef *idef,
113 int ePBC,gmx_bool bMolPBC,t_graph *g,matrix box,
114 t_commrec *cr);
115 /* Spread the force operating on the vsite atoms on the surrounding atoms.
116 * If fshift!=NULL also update the shift forces.
117 * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
118 * to vir. This correction is required when the virial is not calculated
119 * afterwards from the particle position and forces, but in a different way,
120 * as for instance for the PME mesh contribution.
123 GMX_LIBMD_EXPORT
124 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop,t_commrec *cr,
125 gmx_bool bSerial_NoPBC);
126 /* Initialize the virtual site struct,
127 * returns NULL when there are no virtual sites.
128 * bSerial_NoPBC is to generate a simple vsite setup to be
129 * used only serial (no MPI or thread parallelization) and without pbc;
130 * this is useful for correction vsites of the initial configuration.
133 void split_vsites_over_threads(const t_ilist *ilist,
134 const t_mdatoms *mdatoms,
135 gmx_bool bLimitRange,
136 gmx_vsite_t *vsite);
137 /* Divide the vsite work-load over the threads.
138 * Should be called at the end of the domain decomposition.
141 GMX_LIBMD_EXPORT
142 void set_vsite_top(gmx_vsite_t *vsite,gmx_localtop_t *top,t_mdatoms *md,
143 t_commrec *cr);
144 /* Set some vsite data for runs without domain decomposition.
145 * Should be called once after init_vsite, before calling other routines.
148 #ifdef __cplusplus
150 #endif
152 #endif