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