5f93fd7573172eda4259af3bc819fc2cd5c9ec8d
[alexxy/gromacs.git] / src / gromacs / gmxlib / conformation-utilities.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) 2012,2014,2015,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 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "conformation-utilities.h"
41
42 #include <cmath>
43
44 #include <algorithm>
45
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/pbcutil/pbc.h"
49
50 static void low_rotate_conf(int natom, rvec *x, real alfa, real beta, real gamma)
51 {
52     int  i;
53     rvec x_old;
54
55     for (i = 0; i < natom; i++)
56     {
57         copy_rvec(x[i], x_old);
58         /*calculate new x[i] by rotation alfa around the x-axis*/
59         x[i][XX] =   x_old[XX];
60         x[i][YY] =             cos(alfa)*x_old[YY] - sin(alfa)*x_old[ZZ];
61         x[i][ZZ] =             sin(alfa)*x_old[YY] + cos(alfa)*x_old[ZZ];
62         copy_rvec(x[i], x_old);
63         /*calculate new x[i] by rotation beta around the y-axis*/
64         x[i][XX] =   cos(beta)*x_old[XX]           + sin(beta)*x_old[ZZ];
65         x[i][YY] =                       x_old[YY];
66         x[i][ZZ] = -sin(beta)*x_old[XX]           + cos(beta)*x_old[ZZ];
67         copy_rvec(x[i], x_old);
68         /*calculate new x[i] by rotation gamma around the z-axis*/
69         x[i][XX] = x_old[XX]*cos(gamma) - x_old[YY]*sin(gamma);
70         x[i][YY] = x_old[XX]*sin(gamma) + x_old[YY]*cos(gamma);
71         x[i][ZZ] =                                             x_old[ZZ];
72     }
73 }
74
75 void rotate_conf(int natom, rvec *x, rvec *v, real alfa, real beta, real gamma)
76 {
77     if (x)
78     {
79         low_rotate_conf(natom, x, alfa, beta, gamma);
80     }
81     if (v)
82     {
83         low_rotate_conf(natom, v, alfa, beta, gamma);
84     }
85 }
86
87 /* Make a new box around a configuration*/
88 void make_new_box(int natoms, rvec *x, matrix box, const rvec box_space,
89                   gmx_bool bCenter)
90 {
91     int  i, m;
92     rvec xmin, xmax;
93
94     /*calculate minimum and maximum x[0..DIM-1]*/
95     for (m = 0; (m < DIM); m++)
96     {
97         xmin[m] = xmax[m] = x[0][m];
98     }
99     for (i = 1; (i < natoms); i++)
100     {
101         for (m = 0; m < DIM; m++)
102         {
103             xmin[m] = std::min(xmin[m], x[i][m]);
104             xmax[m] = std::max(xmax[m], x[i][m]);
105         }
106     }
107
108     /*calculate the new box sizes for cubic and octahedral ...*/
109     for (m = 0; (m < DIM); m++)
110     {
111         box[m][m] = xmax[m]-xmin[m]+2*box_space[m];
112     }
113
114     /*move the molecule to the center of the box*/
115     if (bCenter)
116     {
117         for (i = 0; (i < natoms); i++)
118         {
119             for (m = 0; (m < DIM); m++)
120             {
121                 x[i][m] += 0.5*(box[m][m]-xmin[m]-xmax[m]);
122             }
123         }
124     }
125 }