Valgrind suppression for OS X 10.9
[alexxy/gromacs.git] / src / programs / view / 3dview.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) 2010,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 "3dview.h"
40
41 #include <math.h>
42
43 #include <algorithm>
44
45 #include "gromacs/math/3dtransforms.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/pbcutil/pbc.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
51
52 static void set_scale(t_3dview *view, real sx, real sy)
53 {
54     view->sc_x = sx;
55     view->sc_y = sy;
56 }
57
58 static void calculate_view(t_3dview *view)
59 {
60 #define SMALL 1e-6
61     mat4 To, Te, T1, T2, T3, T4, T5, N1, D1, D2, D3, D4, D5;
62     real dx, dy, dz, l, r;
63
64     /* eye center */
65     dx = view->eye[XX];
66     dy = view->eye[YY];
67     dz = view->eye[ZZ];
68     l  = sqrt(dx*dx+dy*dy+dz*dz);
69     r  = sqrt(dx*dx+dy*dy);
70 #ifdef DEBUG
71     gmx_vec4_print(debug, "eye", view->eye);
72     printf("del: %10.5f%10.5f%10.5f l: %10.5f, r: %10.5f\n", dx, dy, dz, l, r);
73 #endif
74     if (l < SMALL)
75     {
76         gmx_fatal(FARGS, "Error: Zero Length Vector - No View Specified");
77     }
78     gmx_mat4_init_translation(-view->origin[XX], -view->origin[YY], -view->origin[ZZ], To);
79     gmx_mat4_init_translation(-view->eye[XX], -view->eye[YY], -view->eye[ZZ], Te);
80
81     gmx_mat4_init_unity(T2);
82     T2[YY][YY] = 0, T2[YY][ZZ] = -1, T2[ZZ][YY] = 1, T2[ZZ][ZZ] = 0;
83
84     gmx_mat4_init_unity(T3);
85     if (r > 0)
86     {
87         T3[XX][XX] = -dy/r, T3[XX][ZZ] = dx/r, T3[ZZ][XX] = -dx/r, T3[ZZ][ZZ] = -dy/r;
88     }
89
90     gmx_mat4_init_unity(T4);
91     T4[YY][YY] = r/l, T4[YY][ZZ] = dz/l, T4[ZZ][YY] = -dz/l, T4[ZZ][ZZ] = r/l;
92
93     gmx_mat4_init_unity(T5);
94     T5[ZZ][ZZ] = -1;
95
96     gmx_mat4_init_unity(N1);
97     /* N1[XX][XX]=4,N1[YY][YY]=4; */
98
99     gmx_mat4_mmul(T1, To, view->Rot);
100     gmx_mat4_mmul(D1, Te, T2);
101     gmx_mat4_mmul(D2, T3, T4);
102     gmx_mat4_mmul(D3, T5, N1);
103     gmx_mat4_mmul(D4, T1, D1);
104     gmx_mat4_mmul(D5, D2, D3);
105
106     gmx_mat4_mmul(view->proj, D4, D5);
107
108 #ifdef DEBUG
109     gmx_mat4_print(debug, "T1", T1);
110     gmx_mat4_print(debug, "T2", T2);
111     gmx_mat4_print(debug, "T3", T3);
112     gmx_mat4_print(debug, "T4", T4);
113     gmx_mat4_print(debug, "T5", T5);
114     gmx_mat4_print(debug, "N1", N1);
115     gmx_mat4_print(debug, "Rot", view->Rot);
116     gmx_mat4_print(debug, "Proj", view->proj);
117 #endif
118 }
119
120 gmx_bool zoom_3d(t_3dview *view, real fac)
121 {
122     real dr;
123     real bm, dr1, dr2;
124     int  i;
125
126     dr2 = 0;
127     for (i = 0; (i < DIM); i++)
128     {
129         dr   = view->eye[i];
130         dr2 += dr*dr;
131     }
132     dr1 = sqrt(dr2);
133     if (fac < 1)
134     {
135         bm = std::max(norm(view->box[XX]), std::max(norm(view->box[YY]), norm(view->box[ZZ])));
136         if (dr1*fac < 1.1*bm) /* Don't come to close */
137         {
138             return FALSE;
139         }
140     }
141
142     for (i = 0; (i < DIM); i++)
143     {
144         view->eye[i] *= fac;
145     }
146     calculate_view(view);
147     return TRUE;
148 }
149
150 /* Initiates the state of 3d rotation matrices in the structure */
151 static void init_rotate_3d(t_3dview *view)
152 {
153     real rot = DEG2RAD*15;
154     int  i;
155
156     for (i = 0; (i < DIM); i++)
157     {
158         gmx_mat4_init_rotation(i,  rot, view->RotP[i]);
159         gmx_mat4_init_rotation(i, -rot, view->RotM[i]);
160 #ifdef DEBUG
161         gmx_mat4_print(debug, "RotP", view->RotP[i]);
162         gmx_mat4_print(debug, "RotM", view->RotM[i]);
163 #endif
164     }
165 }
166
167
168 void rotate_3d(t_3dview *view, int axis, gmx_bool bPositive)
169 {
170     mat4 m4;
171
172     if (bPositive)
173     {
174         gmx_mat4_mmul(m4, view->Rot, view->RotP[axis]);
175     }
176     else
177     {
178         gmx_mat4_mmul(m4, view->Rot, view->RotM[axis]);
179     }
180     gmx_mat4_copy(m4, view->Rot);
181     calculate_view(view);
182 }
183
184 void translate_view(t_3dview *view, int axis, gmx_bool bPositive)
185 {
186 #ifdef DEBUG
187     printf("Translate called\n");
188 #endif
189     if (bPositive)
190     {
191         view->origin[axis] += view->box[axis][axis]/8;
192     }
193     else
194     {
195         view->origin[axis] -= view->box[axis][axis]/8;
196     }
197     calculate_view(view);
198 }
199
200 void reset_view(t_3dview *view)
201 {
202 #ifdef DEBUG
203     printf("Reset view called\n");
204 #endif
205     set_scale(view, 4.0, 4.0);
206     clear_rvec(view->eye);
207     calc_box_center(view->ecenter, view->box, view->origin);
208     view->eye[ZZ] = 3.0*std::max(view->box[XX][XX], view->box[YY][YY]);
209     zoom_3d(view, 1.0);
210     view->eye[WW] = view->origin[WW] = 0.0;
211
212     /* Initiate the matrix */
213     gmx_mat4_init_unity(view->Rot);
214     calculate_view(view);
215
216     init_rotate_3d(view);
217 }
218
219 t_3dview *init_view(matrix box)
220 {
221     t_3dview *view;
222
223     snew(view, 1);
224     copy_mat(box, view->box);
225     view->ecenter = ecenterDEF;
226     reset_view(view);
227
228     return view;
229 }