Remove continuation from convert_tpr
[gromacs.git] / src / gromacs / pbcutil / mshift.h
blobc5540bc34b8647b20106b9793e5be8c37cdcd03f
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_PBCUTIL_MSHIFT_H
38 #define GMX_PBCUTIL_MSHIFT_H
40 #include <stdio.h>
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/basedefinitions.h"
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
49 struct t_idef;
50 struct t_ilist;
52 typedef enum {
53 egcolWhite, egcolGrey, egcolBlack, egcolNR
54 } egCol;
56 typedef struct t_graph {
57 int at0; /* The first atom the graph was constructed for */
58 int at1; /* The last atom the graph was constructed for */
59 int nnodes; /* The number of nodes, nnodes=at_end-at_start */
60 int nbound; /* The number of nodes with edges */
61 int at_start; /* The first connected atom in this graph */
62 int at_end; /* The last+1 connected atom in this graph */
63 int *nedge; /* For each node the number of edges */
64 int **edge; /* For each node, the actual edges (bidirect.) */
65 gmx_bool bScrewPBC; /* Screw boundary conditions */
66 ivec *ishift; /* Shift for each particle */
67 int negc;
68 egCol *egc; /* color of each node */
69 } t_graph;
71 #define SHIFT_IVEC(g, i) ((g)->ishift[i])
73 t_graph *mk_graph(FILE *fplog,
74 const struct t_idef *idef, int at_start, int at_end,
75 gmx_bool bShakeOnly, gmx_bool bSettle);
76 /* Build a graph from an idef description. The graph can be used
77 * to generate mol-shift indices.
78 * at_start and at_end should coincide will molecule boundaries,
79 * for the whole system this is simply 0 and natoms.
80 * If bShakeOnly, only the connections in the shake list are used.
81 * If bSettle && bShakeOnly the settles are used too.
84 void mk_graph_ilist(FILE *fplog,
85 const struct t_ilist *ilist, int at_start, int at_end,
86 gmx_bool bShakeOnly, gmx_bool bSettle,
87 t_graph *g);
88 /* As mk_graph, but takes t_ilist iso t_idef and does not allocate g */
91 void done_graph(t_graph *g);
92 /* Free the memory in g */
94 void p_graph(FILE *log, const char *title, t_graph *g);
95 /* Print a graph to log */
97 void mk_mshift(FILE *log, t_graph *g, int ePBC,
98 const matrix box, const rvec x[]);
99 /* Calculate the mshift codes, based on the connection graph in g. */
101 void shift_x(const t_graph *g, const matrix box, const rvec x[], rvec x_s[]);
102 /* Add the shift vector to x, and store in x_s (may be same array as x) */
104 void shift_self(const t_graph *g, const matrix box, rvec x[]);
105 /* Id. but in place */
107 void unshift_x(const t_graph *g, const matrix box, rvec x[], const rvec x_s[]);
108 /* Subtract the shift vector from x_s, and store in x (may be same array) */
110 void unshift_self(const t_graph *g, const matrix box, rvec x[]);
111 /* Id, but in place */
113 #ifdef __cplusplus
115 #endif
117 #endif