Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxlib / conformation-utilities.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) 2012,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 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "conformation-utilities.h"
41
42 #include <math.h>
43
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/pbcutil/pbc.h"
47
48 static real dist2(t_pbc *pbc, rvec x, rvec y)
49 {
50     rvec dx;
51
52     pbc_dx(pbc, x, y, dx);
53
54     return norm2(dx);
55 }
56
57 static real distance_to_z(rvec x)
58 {
59     return (sqr(x[XX])+sqr(x[YY]));
60 } /*distance_to_z()*/
61
62 static void low_rotate_conf(int natom, rvec *x, real alfa, real beta, real gamma)
63 {
64     int  i;
65     rvec x_old;
66
67     for (i = 0; i < natom; i++)
68     {
69         copy_rvec(x[i], x_old);
70         /*calculate new x[i] by rotation alfa around the x-axis*/
71         x[i][XX] =   x_old[XX];
72         x[i][YY] =             cos(alfa)*x_old[YY] - sin(alfa)*x_old[ZZ];
73         x[i][ZZ] =             sin(alfa)*x_old[YY] + cos(alfa)*x_old[ZZ];
74         copy_rvec(x[i], x_old);
75         /*calculate new x[i] by rotation beta around the y-axis*/
76         x[i][XX] =   cos(beta)*x_old[XX]           + sin(beta)*x_old[ZZ];
77         x[i][YY] =                       x_old[YY];
78         x[i][ZZ] = -sin(beta)*x_old[XX]           + cos(beta)*x_old[ZZ];
79         copy_rvec(x[i], x_old);
80         /*calculate new x[i] by rotation gamma around the z-axis*/
81         x[i][XX] = x_old[XX]*cos(gamma) - x_old[YY]*sin(gamma);
82         x[i][YY] = x_old[XX]*sin(gamma) + x_old[YY]*cos(gamma);
83         x[i][ZZ] =                                             x_old[ZZ];
84     }
85 }
86
87 static void low_rotate_conf_indexed(int nindex, atom_id *index, rvec *x, real alfa, real beta, real gamma)
88 {
89     int  i;
90     rvec x_old;
91
92     for (i = 0; i < nindex; i++)
93     {
94         copy_rvec(x[index[i]], x_old);
95         /*calculate new x[index[i]] by rotation alfa around the x-axis*/
96         x[index[i]][XX] =   x_old[XX];
97         x[index[i]][YY] =             cos(alfa)*x_old[YY] - sin(alfa)*x_old[ZZ];
98         x[index[i]][ZZ] =             sin(alfa)*x_old[YY] + cos(alfa)*x_old[ZZ];
99         copy_rvec(x[index[i]], x_old);
100         /*calculate new x[index[i]] by rotation beta around the y-axis*/
101         x[index[i]][XX] =   cos(beta)*x_old[XX]           + sin(beta)*x_old[ZZ];
102         x[index[i]][YY] =                       x_old[YY];
103         x[index[i]][ZZ] = -sin(beta)*x_old[XX]           + cos(beta)*x_old[ZZ];
104         copy_rvec(x[index[i]], x_old);
105         /*calculate new x[index[i]] by rotation gamma around the z-axis*/
106         x[index[i]][XX] = x_old[XX]*cos(gamma) - x_old[YY]*sin(gamma);
107         x[index[i]][YY] = x_old[XX]*sin(gamma) + x_old[YY]*cos(gamma);
108         x[index[i]][ZZ] =                                             x_old[ZZ];
109     }
110 }
111
112 void rotate_conf(int natom, rvec *x, rvec *v, real alfa, real beta, real gamma)
113 {
114     if (x)
115     {
116         low_rotate_conf(natom, x, alfa, beta, gamma);
117     }
118     if (v)
119     {
120         low_rotate_conf(natom, v, alfa, beta, gamma);
121     }
122 }
123
124 /* Make a new box around a configuration*/
125 void make_new_box(int natoms, rvec *x, matrix box, rvec box_space,
126                   gmx_bool bCenter)
127 {
128     int  i, m;
129     rvec xmin, xmax;
130     real max_box;
131
132     /*calculate minimum and maximum x[0..DIM-1]*/
133     for (m = 0; (m < DIM); m++)
134     {
135         xmin[m] = xmax[m] = x[0][m];
136     }
137     for (i = 1; (i < natoms); i++)
138     {
139         for (m = 0; m < DIM; m++)
140         {
141             xmin[m] = min(xmin[m], x[i][m]);
142             xmax[m] = max(xmax[m], x[i][m]);
143         }
144     }
145
146     /*calculate the new box sizes for cubic and octahedral ...*/
147     for (m = 0; (m < DIM); m++)
148     {
149         box[m][m] = xmax[m]-xmin[m]+2*box_space[m];
150     }
151
152     /*move the molecule to the center of the box*/
153     if (bCenter)
154     {
155         for (i = 0; (i < natoms); i++)
156         {
157             for (m = 0; (m < DIM); m++)
158             {
159                 x[i][m] += 0.5*(box[m][m]-xmin[m]-xmax[m]);
160             }
161         }
162     }
163 }