Only issue FFT warning messages on changes
[gromacs/AngularHB.git] / src / programs / view / 3dview.cpp
blobaf36b05cad1bbb7f9d02278e42536b595055f80c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2010,2014,2015, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "3dview.h"
41 #include <cmath>
42 #include <cstdio>
44 #include <algorithm>
46 #include "gromacs/math/3dtransforms.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 static void set_scale(t_3dview *view, real sx, real sy)
55 view->sc_x = sx;
56 view->sc_y = sy;
59 static void calculate_view(t_3dview *view)
61 #define SMALL 1e-6
62 mat4 To, Te, T1, T2, T3, T4, T5, N1, D1, D2, D3, D4, D5;
63 real dx, dy, dz, l, r;
65 /* eye center */
66 dx = view->eye[XX];
67 dy = view->eye[YY];
68 dz = view->eye[ZZ];
69 l = std::sqrt(dx*dx+dy*dy+dz*dz);
70 r = std::sqrt(dx*dx+dy*dy);
71 #ifdef DEBUG
72 gmx_vec4_print(debug, "eye", view->eye);
73 std::printf("del: %10.5f%10.5f%10.5f l: %10.5f, r: %10.5f\n", dx, dy, dz, l, r);
74 #endif
75 if (l < SMALL)
77 gmx_fatal(FARGS, "Error: Zero Length Vector - No View Specified");
79 gmx_mat4_init_translation(-view->origin[XX], -view->origin[YY], -view->origin[ZZ], To);
80 gmx_mat4_init_translation(-view->eye[XX], -view->eye[YY], -view->eye[ZZ], Te);
82 gmx_mat4_init_unity(T2);
83 T2[YY][YY] = 0, T2[YY][ZZ] = -1, T2[ZZ][YY] = 1, T2[ZZ][ZZ] = 0;
85 gmx_mat4_init_unity(T3);
86 if (r > 0)
88 T3[XX][XX] = -dy/r, T3[XX][ZZ] = dx/r, T3[ZZ][XX] = -dx/r, T3[ZZ][ZZ] = -dy/r;
91 gmx_mat4_init_unity(T4);
92 T4[YY][YY] = r/l, T4[YY][ZZ] = dz/l, T4[ZZ][YY] = -dz/l, T4[ZZ][ZZ] = r/l;
94 gmx_mat4_init_unity(T5);
95 T5[ZZ][ZZ] = -1;
97 gmx_mat4_init_unity(N1);
98 /* N1[XX][XX]=4,N1[YY][YY]=4; */
100 gmx_mat4_mmul(T1, To, view->Rot);
101 gmx_mat4_mmul(D1, Te, T2);
102 gmx_mat4_mmul(D2, T3, T4);
103 gmx_mat4_mmul(D3, T5, N1);
104 gmx_mat4_mmul(D4, T1, D1);
105 gmx_mat4_mmul(D5, D2, D3);
107 gmx_mat4_mmul(view->proj, D4, D5);
109 #ifdef DEBUG
110 gmx_mat4_print(debug, "T1", T1);
111 gmx_mat4_print(debug, "T2", T2);
112 gmx_mat4_print(debug, "T3", T3);
113 gmx_mat4_print(debug, "T4", T4);
114 gmx_mat4_print(debug, "T5", T5);
115 gmx_mat4_print(debug, "N1", N1);
116 gmx_mat4_print(debug, "Rot", view->Rot);
117 gmx_mat4_print(debug, "Proj", view->proj);
118 #endif
121 gmx_bool zoom_3d(t_3dview *view, real fac)
123 real dr;
124 real bm, dr1, dr2;
125 int i;
127 dr2 = 0;
128 for (i = 0; (i < DIM); i++)
130 dr = view->eye[i];
131 dr2 += dr*dr;
133 dr1 = std::sqrt(dr2);
134 if (fac < 1)
136 bm = std::max(norm(view->box[XX]), std::max(norm(view->box[YY]), norm(view->box[ZZ])));
137 if (dr1*fac < 1.1*bm) /* Don't come to close */
139 return FALSE;
143 for (i = 0; (i < DIM); i++)
145 view->eye[i] *= fac;
147 calculate_view(view);
148 return TRUE;
151 /* Initiates the state of 3d rotation matrices in the structure */
152 static void init_rotate_3d(t_3dview *view)
154 real rot = DEG2RAD*15;
155 int i;
157 for (i = 0; (i < DIM); i++)
159 gmx_mat4_init_rotation(i, rot, view->RotP[i]);
160 gmx_mat4_init_rotation(i, -rot, view->RotM[i]);
161 #ifdef DEBUG
162 gmx_mat4_print(debug, "RotP", view->RotP[i]);
163 gmx_mat4_print(debug, "RotM", view->RotM[i]);
164 #endif
169 void rotate_3d(t_3dview *view, int axis, gmx_bool bPositive)
171 mat4 m4;
173 if (bPositive)
175 gmx_mat4_mmul(m4, view->Rot, view->RotP[axis]);
177 else
179 gmx_mat4_mmul(m4, view->Rot, view->RotM[axis]);
181 gmx_mat4_copy(m4, view->Rot);
182 calculate_view(view);
185 void translate_view(t_3dview *view, int axis, gmx_bool bPositive)
187 #ifdef DEBUG
188 std::printf("Translate called\n");
189 #endif
190 if (bPositive)
192 view->origin[axis] += view->box[axis][axis]/8;
194 else
196 view->origin[axis] -= view->box[axis][axis]/8;
198 calculate_view(view);
201 void reset_view(t_3dview *view)
203 #ifdef DEBUG
204 std::printf("Reset view called\n");
205 #endif
206 set_scale(view, 4.0, 4.0);
207 clear_rvec(view->eye);
208 calc_box_center(view->ecenter, view->box, view->origin);
209 view->eye[ZZ] = 3.0*std::max(view->box[XX][XX], view->box[YY][YY]);
210 zoom_3d(view, 1.0);
211 view->eye[WW] = view->origin[WW] = 0.0;
213 /* Initiate the matrix */
214 gmx_mat4_init_unity(view->Rot);
215 calculate_view(view);
217 init_rotate_3d(view);
220 t_3dview *init_view(matrix box)
222 t_3dview *view;
224 snew(view, 1);
225 copy_mat(box, view->box);
226 view->ecenter = ecenterDEF;
227 reset_view(view);
229 return view;