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