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