added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / include / vsite.h
blob67dc8207010c7c697a0829e8d3a3a64a039742db
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
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
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _vsite_h
37 #define _vsite_h
39 #include <stdio.h>
40 #include "typedefs.h"
41 #include "types/commrec.h"
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
47 typedef struct {
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;
56 rvec * send_buf;
57 rvec * recv_buf;
58 } t_comm_vsites;
60 typedef struct {
61 int n_intercg_vsite; /* The number of inter charge group vsites */
62 int nvsite_pbc_molt; /* The array size of vsite_pbc_molt */
63 int ***vsite_pbc_molt; /* The pbc atoms for intercg vsites */
64 int **vsite_pbc_loc; /* The local pbc atoms */
65 int *vsite_pbc_loc_nalloc;
66 gmx_bool bPDvsitecomm; /* Do we need vsite communication with PD? */
67 t_comm_vsites *vsitecomm; /* The PD vsite communication struct */
68 } gmx_vsite_t;
70 void construct_vsites(FILE *log,gmx_vsite_t *vsite,
71 rvec x[],t_nrnb *nrnb,
72 real dt,rvec v[],
73 t_iparams ip[],t_ilist ilist[],
74 int ePBC,gmx_bool bMolPBC,t_graph *graph,
75 t_commrec *cr,matrix box);
76 /* Create positions of vsite atoms based on surrounding atoms
77 * for the local system.
78 * If v is passed, the velocities of the vsites will be calculated
79 * as the new positions minus the old positions divided by dt,
80 * thus v should only be passed when the coordinates have been
81 * updated with a full time step.
82 * Note that velocitis of vsites are completely irrelevant
83 * for the integration, they are only useful for analysis.
86 void construct_vsites_mtop(FILE *log,gmx_vsite_t *vsite,
87 gmx_mtop_t *mtop,rvec x[]);
88 /* Create positions of vsite atoms based on surrounding atoms
89 * for the whole system.
90 * This function assumes that all molecules are whole.
93 void spread_vsite_f(FILE *log,gmx_vsite_t *vsite,
94 rvec x[],rvec f[],rvec *fshift,
95 gmx_bool VirCorr,matrix vir,
96 t_nrnb *nrnb,t_idef *idef,
97 int ePBC,gmx_bool bMolPBC,t_graph *g,matrix box,
98 t_commrec *cr);
99 /* Spread the force operating on the vsite atoms on the surrounding atoms.
100 * If fshift!=NULL also update the shift forces.
101 * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
102 * to vir. This correction is required when the virial is not calculated
103 * afterwards from the particle position and forces, but in a different way,
104 * as for instance for the PME mesh contribution.
107 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop,t_commrec *cr);
108 /* Initialize the virtual site struct,
109 * returns NULL when there are no virtual sites.
112 void set_vsite_top(gmx_vsite_t *vsite,gmx_localtop_t *top,t_mdatoms *md,
113 t_commrec *cr);
114 /* Set some vsite data for runs without domain decomposition.
115 * Should be called once after init_vsite, before calling other routines.
118 #ifdef __cplusplus
120 #endif
122 #endif