Merge branch release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxlib / rmpbc.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "sysstuff.h"
42 #include "typedefs.h"
43 #include "mshift.h"
44 #include "pbc.h"
45 #include "rmpbc.h"
46 #include "vec.h"
47
48 #include "gromacs/fileio/futil.h"
49 #include "gmx_fatal.h"
50 #include "gromacs/utility/smalloc.h"
51
52 typedef struct {
53     int      natoms;
54     t_graph *gr;
55 } rmpbc_graph_t;
56
57 typedef struct gmx_rmpbc {
58     t_idef        *idef;
59     int            natoms_init;
60     int            ePBC;
61     int            ngraph;
62     rmpbc_graph_t *graph;
63 } koeiepoep;
64
65 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc, int ePBC, int natoms)
66 {
67     int            i;
68     rmpbc_graph_t *gr;
69
70     if (ePBC == epbcNONE
71         || NULL == gpbc
72         || NULL == gpbc->idef
73         || gpbc->idef->ntypes <= 0)
74     {
75         return NULL;
76     }
77
78     gr = NULL;
79     for (i = 0; i < gpbc->ngraph; i++)
80     {
81         if (natoms == gpbc->graph[i].natoms)
82         {
83             gr = &gpbc->graph[i];
84         }
85     }
86     if (gr == NULL)
87     {
88         /* We'd like to check with the number of atoms in the topology,
89          * but we don't have that available.
90          * So we check against the number of atoms that gmx_rmpbc_init
91          * was called with.
92          */
93         if (natoms > gpbc->natoms_init)
94         {
95             gmx_fatal(FARGS, "Structure or trajectory file has more atoms (%d) than the topology (%d)", natoms, gpbc->natoms_init);
96         }
97         gpbc->ngraph++;
98         srenew(gpbc->graph, gpbc->ngraph);
99         gr         = &gpbc->graph[gpbc->ngraph-1];
100         gr->natoms = natoms;
101         gr->gr     = mk_graph(NULL, gpbc->idef, 0, natoms, FALSE, FALSE);
102     }
103
104     return gr->gr;
105 }
106
107 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef, int ePBC, int natoms)
108 {
109     gmx_rmpbc_t gpbc;
110
111     snew(gpbc, 1);
112
113     gpbc->natoms_init = natoms;
114
115     /* This sets pbc when we now it,
116      * otherwise we guess it from the instantaneous box in the trajectory.
117      */
118     gpbc->ePBC = ePBC;
119
120     gpbc->idef = idef;
121     if (gpbc->idef->ntypes <= 0)
122     {
123         fprintf(stderr,
124                 "\n"
125                 "WARNING: If there are molecules in the input trajectory file\n"
126                 "         that are broken across periodic boundaries, they\n"
127                 "         cannot be made whole (or treated as whole) without\n"
128                 "         you providing a run input file.\n\n");
129     }
130
131     return gpbc;
132 }
133
134 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
135 {
136     int i;
137
138     if (NULL != gpbc)
139     {
140         for (i = 0; i < gpbc->ngraph; i++)
141         {
142             done_graph(gpbc->graph[i].gr);
143             sfree(gpbc->graph[i].gr);
144         }
145         if (gpbc->graph != NULL)
146         {
147             sfree(gpbc->graph);
148         }
149         sfree(gpbc);
150     }
151 }
152
153 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc, matrix box)
154 {
155     if (NULL != gpbc && gpbc->ePBC >= 0)
156     {
157         return gpbc->ePBC;
158     }
159     else
160     {
161         return guess_ePBC(box);
162     }
163 }
164
165 void gmx_rmpbc(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[])
166 {
167     int      ePBC;
168     t_graph *gr;
169
170     ePBC = gmx_rmpbc_ePBC(gpbc, box);
171     gr   = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
172     if (gr != NULL)
173     {
174         mk_mshift(stdout, gr, ePBC, box, x);
175         shift_self(gr, box, x);
176     }
177 }
178
179 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[], rvec x_s[])
180 {
181     int      ePBC;
182     t_graph *gr;
183     int      i;
184
185     ePBC = gmx_rmpbc_ePBC(gpbc, box);
186     gr   = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
187     if (gr != NULL)
188     {
189         mk_mshift(stdout, gr, ePBC, box, x);
190         shift_x(gr, box, x, x_s);
191     }
192     else
193     {
194         for (i = 0; i < natoms; i++)
195         {
196             copy_rvec(x[i], x_s[i]);
197         }
198     }
199 }
200
201 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc, t_trxframe *fr)
202 {
203     int      ePBC;
204     t_graph *gr;
205
206     if (fr->bX && fr->bBox)
207     {
208         ePBC = gmx_rmpbc_ePBC(gpbc, fr->box);
209         gr   = gmx_rmpbc_get_graph(gpbc, ePBC, fr->natoms);
210         if (gr != NULL)
211         {
212             mk_mshift(stdout, gr, ePBC, fr->box, fr->x);
213             shift_self(gr, fr->box, fr->x);
214         }
215     }
216 }
217
218 void rm_gropbc(t_atoms *atoms, rvec x[], matrix box)
219 {
220     real dist;
221     int  n, m, d;
222
223     /* check periodic boundary */
224     for (n = 1; (n < atoms->nr); n++)
225     {
226         for (m = DIM-1; m >= 0; m--)
227         {
228             dist = x[n][m]-x[n-1][m];
229             if (fabs(dist) > 0.9*box[m][m])
230             {
231                 if (dist >  0)
232                 {
233                     for (d = 0; d <= m; d++)
234                     {
235                         x[n][d] -= box[m][d];
236                     }
237                 }
238                 else
239                 {
240                     for (d = 0; d <= m; d++)
241                     {
242                         x[n][d] += box[m][d];
243                     }
244                 }
245             }
246         }
247     }
248 }