70b70dc632d109103707fb08be31438f385b3698
[alexxy/gromacs.git] / src / gromacs / gmxlib / rmpbc.c
1 /*  -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "sysstuff.h"
41 #include "typedefs.h"
42 #include "smalloc.h"
43 #include "mshift.h"
44 #include "pbc.h"
45 #include "gstat.h"
46 #include "futil.h"
47 #include "vec.h"
48
49 typedef struct {
50     int      natoms;
51     t_graph *gr;
52 } rmpbc_graph_t;
53
54 typedef struct gmx_rmpbc {
55     t_idef        *idef;
56     int            natoms_init;
57     int            ePBC;
58     int            ngraph;
59     rmpbc_graph_t *graph;
60 } koeiepoep;
61
62 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc, int ePBC, int natoms)
63 {
64     int            i;
65     rmpbc_graph_t *gr;
66
67     if (ePBC == epbcNONE
68         || NULL == gpbc
69         || NULL == gpbc->idef
70         || gpbc->idef->ntypes <= 0)
71     {
72         return NULL;
73     }
74
75     gr = NULL;
76     for (i = 0; i < gpbc->ngraph; i++)
77     {
78         if (natoms == gpbc->graph[i].natoms)
79         {
80             gr = &gpbc->graph[i];
81         }
82     }
83     if (gr == NULL)
84     {
85         /* We'd like to check with the number of atoms in the topology,
86          * but we don't have that available.
87          * So we check against the number of atoms that gmx_rmpbc_init
88          * was called with.
89          */
90         if (natoms > gpbc->natoms_init)
91         {
92             gmx_fatal(FARGS, "Structure or trajectory file has more atoms (%d) than the topology (%d)", natoms, gpbc->natoms_init);
93         }
94         gpbc->ngraph++;
95         srenew(gpbc->graph, gpbc->ngraph);
96         gr         = &gpbc->graph[gpbc->ngraph-1];
97         gr->natoms = natoms;
98         gr->gr     = mk_graph(NULL, gpbc->idef, 0, natoms, FALSE, FALSE);
99     }
100
101     return gr->gr;
102 }
103
104 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef, int ePBC, int natoms,
105                            matrix box)
106 {
107     gmx_rmpbc_t gpbc;
108
109     snew(gpbc, 1);
110
111     gpbc->natoms_init = natoms;
112
113     /* This sets pbc when we now it,
114      * otherwise we guess it from the instantaneous box in the trajectory.
115      */
116     gpbc->ePBC = ePBC;
117
118     gpbc->idef = idef;
119     if (gpbc->idef->ntypes <= 0)
120     {
121         fprintf(stderr,
122                 "\n"
123                 "WARNING: If there are molecules in the input trajectory file\n"
124                 "         that are broken across periodic boundaries, they\n"
125                 "         cannot be made whole (or treated as whole) without\n"
126                 "         you providing a run input file.\n\n");
127     }
128
129     return gpbc;
130 }
131
132 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
133 {
134     int i;
135
136     if (NULL != gpbc)
137     {
138         for (i = 0; i < gpbc->ngraph; i++)
139         {
140             done_graph(gpbc->graph[i].gr);
141         }
142         if (gpbc->graph != NULL)
143         {
144             sfree(gpbc->graph);
145         }
146         sfree(gpbc);
147     }
148 }
149
150 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc, matrix box)
151 {
152     if (NULL != gpbc && gpbc->ePBC >= 0)
153     {
154         return gpbc->ePBC;
155     }
156     else
157     {
158         return guess_ePBC(box);
159     }
160 }
161
162 void gmx_rmpbc(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[])
163 {
164     int      ePBC;
165     t_graph *gr;
166
167     ePBC = gmx_rmpbc_ePBC(gpbc, box);
168     gr   = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
169     if (gr != NULL)
170     {
171         mk_mshift(stdout, gr, ePBC, box, x);
172         shift_self(gr, box, x);
173     }
174 }
175
176 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[], rvec x_s[])
177 {
178     int      ePBC;
179     t_graph *gr;
180     int      i;
181
182     ePBC = gmx_rmpbc_ePBC(gpbc, box);
183     gr   = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
184     if (gr != NULL)
185     {
186         mk_mshift(stdout, gr, ePBC, box, x);
187         shift_x(gr, box, x, x_s);
188     }
189     else
190     {
191         for (i = 0; i < natoms; i++)
192         {
193             copy_rvec(x[i], x_s[i]);
194         }
195     }
196 }
197
198 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc, t_trxframe *fr)
199 {
200     int      ePBC;
201     t_graph *gr;
202
203     if (fr->bX && fr->bBox)
204     {
205         ePBC = gmx_rmpbc_ePBC(gpbc, fr->box);
206         gr   = gmx_rmpbc_get_graph(gpbc, ePBC, fr->natoms);
207         if (gr != NULL)
208         {
209             mk_mshift(stdout, gr, ePBC, fr->box, fr->x);
210             shift_self(gr, fr->box, fr->x);
211         }
212     }
213 }
214
215 void rm_gropbc(t_atoms *atoms, rvec x[], matrix box)
216 {
217     real dist;
218     int  n, m, d;
219
220     /* check periodic boundary */
221     for (n = 1; (n < atoms->nr); n++)
222     {
223         for (m = DIM-1; m >= 0; m--)
224         {
225             dist = x[n][m]-x[n-1][m];
226             if (fabs(dist) > 0.9*box[m][m])
227             {
228                 if (dist >  0)
229                 {
230                     for (d = 0; d <= m; d++)
231                     {
232                         x[n][d] -= box[m][d];
233                     }
234                 }
235                 else
236                 {
237                     for (d = 0; d <= m; d++)
238                     {
239                         x[n][d] += box[m][d];
240                     }
241                 }
242             }
243         }
244     }
245 }