Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / vsite.h
blob051f9b7c31432d8a4c282d5e04cf15551320755d
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 * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MDLIB_VSITE_H
38 #define GMX_MDLIB_VSITE_H
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/pbcutil/ishift.h"
42 #include "gromacs/topology/idef.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
46 struct gmx_localtop_t;
47 struct gmx_mtop_t;
48 struct t_commrec;
49 struct t_graph;
50 struct t_ilist;
51 struct t_mdatoms;
52 struct t_nrnb;
54 typedef struct gmx_vsite_t {
55 gmx_bool bHaveChargeGroups; /* Do we have charge groups? */
56 int n_intercg_vsite; /* The number of inter charge group vsites */
57 int nvsite_pbc_molt; /* The array size of vsite_pbc_molt */
58 int ***vsite_pbc_molt; /* The pbc atoms for intercg vsites */
59 int **vsite_pbc_loc; /* The local pbc atoms */
60 int *vsite_pbc_loc_nalloc; /* Sizes of vsite_pbc_loc */
61 int nthreads; /* Number of threads used for vsites */
62 struct VsiteThread **tData; /* Thread local vsites and work structs */
63 int *taskIndex; /* Work array */
64 int taskIndexNalloc; /* Size of taskIndex */
65 } gmx_vsite_t;
67 void construct_vsites(const gmx_vsite_t *vsite,
68 rvec x[],
69 real dt, rvec v[],
70 const t_iparams ip[], const t_ilist ilist[],
71 int ePBC, gmx_bool bMolPBC,
72 t_commrec *cr, matrix box);
73 /* Create positions of vsite atoms based on surrounding atoms
74 * for the local system.
75 * If v is passed, the velocities of the vsites will be calculated
76 * as the new positions minus the old positions divided by dt,
77 * thus v should only be passed when the coordinates have been
78 * updated with a full time step.
79 * Note that velocitis of vsites are completely irrelevant
80 * for the integration, they are only useful for analysis.
83 void construct_vsites_mtop(gmx_vsite_t *vsite,
84 gmx_mtop_t *mtop, rvec x[]);
85 /* Create positions of vsite atoms based on surrounding atoms
86 * for the whole system.
87 * This function assumes that all molecules are whole.
90 void spread_vsite_f(const gmx_vsite_t *vsite,
91 const rvec x[],
92 rvec f[], rvec *fshift,
93 gmx_bool VirCorr, matrix vir,
94 t_nrnb *nrnb, const t_idef *idef,
95 int ePBC, gmx_bool bMolPBC, const t_graph *g, const matrix box,
96 t_commrec *cr);
97 /* Spread the force operating on the vsite atoms on the surrounding atoms.
98 * If fshift!=NULL also update the shift forces.
99 * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
100 * to vir. This correction is required when the virial is not calculated
101 * afterwards from the particle position and forces, but in a different way,
102 * as for instance for the PME mesh contribution.
105 int count_intercg_vsites(const gmx_mtop_t *mtop);
106 /* Returns the number of virtual sites that cross charge groups */
108 gmx_vsite_t *init_vsite(const gmx_mtop_t *mtop, t_commrec *cr,
109 gmx_bool bSerial_NoPBC);
110 /* Initialize the virtual site struct,
111 * returns NULL when there are no virtual sites.
112 * bSerial_NoPBC is to generate a simple vsite setup to be
113 * used only serial (no MPI or thread parallelization) and without pbc;
114 * this is useful for correction vsites of the initial configuration.
117 void split_vsites_over_threads(const t_ilist *ilist,
118 const t_iparams *ip,
119 const t_mdatoms *mdatoms,
120 gmx_bool bLimitRange,
121 gmx_vsite_t *vsite);
122 /* Divide the vsite work-load over the threads.
123 * Should be called at the end of the domain decomposition.
126 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
127 t_commrec *cr);
128 /* Set some vsite data for runs without domain decomposition.
129 * Should be called once after init_vsite, before calling other routines.
132 #endif