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