Split txtdump.*, part 2
[gromacs/AngularHB.git] / src / programs / view / nmol.cpp
blobd2f27ef9d3f9b609d1853058ed92ebcd5d1e3488
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-2013, The GROMACS development team.
6 * Copyright (c) 2013,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 "nmol.h"
41 #include <cstdio>
42 #include <cstdlib>
43 #include <cstring>
45 #include "gromacs/math/vec.h"
46 #include "gromacs/math/vecdump.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
51 #include "3dview.h"
52 #include "buttons.h"
53 #include "manager.h"
54 #include "xutil.h"
56 #define MSIZE 4
58 static bool MWCallBack(t_x11 *x11, XEvent *event, Window /*w*/, void *data)
60 t_molwin *mw;
61 Window To;
62 XEvent letter;
64 mw = (t_molwin *)data;
65 To = mw->wd.Parent;
66 letter.type = ClientMessage;
67 letter.xclient.display = x11->disp;
68 letter.xclient.window = To;
69 letter.xclient.message_type = 0;
70 letter.xclient.format = 32;
71 switch (event->type)
73 case Expose:
74 /* Do Not draw anything, but signal parent instead, he will
75 * coordinate drawing.
77 letter.xclient.data.l[0] = IDDRAWMOL;
78 letter.xclient.data.l[1] = Button1;
79 XSendEvent(x11->disp, To, True, 0, &letter);
80 break;
81 case ButtonPress:
82 #ifdef DEBUG
83 std::printf("Molwindow: Buttonpress\n");
84 #endif
85 letter.xclient.data.l[0] = IDLABEL;
86 letter.xclient.data.l[1] = (long)event->xbutton.button;
87 letter.xclient.data.l[2] = event->xbutton.x;
88 letter.xclient.data.l[3] = event->xbutton.y;
89 XSendEvent(x11->disp, To, True, 0, &letter);
90 break;
91 case ConfigureNotify:
92 mw->wd.width = event->xconfigure.width;
93 mw->wd.height = event->xconfigure.height;
94 break;
95 default:
96 break;
98 return false;
101 void set_def (t_molwin *mw, int ePBC, matrix box)
103 mw->bShowHydrogen = true;
104 mw->bond_type = eBFat;
105 mw->ePBC = ePBC;
106 mw->boxtype = esbRect;
107 mw->realbox = TRICLINIC(box) ? esbTri : esbRect;
110 t_molwin *init_mw(t_x11 *x11, Window Parent,
111 int x, int y, int width, int height,
112 unsigned long fg, unsigned long bg,
113 int ePBC, matrix box)
115 t_molwin *mw;
117 snew(mw, 1);
118 set_def(mw, ePBC, box);
120 InitWin(&mw->wd, x, y, width, height, 1, "Mol Window");
122 mw->wd.Parent = Parent;
123 mw->wd.self = XCreateSimpleWindow(x11->disp, Parent, x, y, width, height, 1, fg, bg);
124 x11->RegisterCallback(x11, mw->wd.self, Parent, MWCallBack, mw);
125 x11->SetInputMask(x11, mw->wd.self,
126 ExposureMask | StructureNotifyMask |
127 ButtonPressMask);
128 return mw;
131 void map_mw(t_x11 *x11, t_molwin *mw)
133 XMapWindow(x11->disp, mw->wd.self);
136 bool toggle_hydrogen(t_x11 *x11, t_molwin *mw)
138 mw->bShowHydrogen = !mw->bShowHydrogen;
139 ExposeWin(x11->disp, mw->wd.self);
141 return mw->bShowHydrogen;
144 void set_bond_type(t_x11 *x11, t_molwin *mw, int bt)
146 if (bt != mw->bond_type)
148 mw->bond_type = bt;
149 ExposeWin(x11->disp, mw->wd.self);
153 void set_box_type (t_x11 *x11, t_molwin *mw, int bt)
155 #ifdef DEBUG
156 std::fprintf(stderr, "mw->boxtype = %d, bt = %d\n", mw->boxtype, bt);
157 #endif
158 if (bt != mw->boxtype)
160 if ((bt == esbTrunc && mw->realbox == esbTri) || bt == esbTri || bt == esbNone)
162 mw->boxtype = bt;
163 ExposeWin(x11->disp, mw->wd.self);
165 else
167 std::fprintf(stderr, "Can not change rectangular box to truncated octahedron\n");
172 void done_mw(t_x11 *x11, t_molwin *mw)
174 x11->UnRegisterCallback(x11, mw->wd.self);
175 sfree(mw);
178 static void draw_atom(Display *disp, Window w, GC gc,
179 int ai, iv2 vec2[], unsigned long col[], int size[],
180 bool bBall, bool bPlus)
182 int xi, yi;
184 xi = vec2[ai][XX];
185 yi = vec2[ai][YY];
186 XSetForeground(disp, gc, col[ai]);
187 if (bBall)
189 XFillCircle(disp, w, gc, xi, yi, size[ai]-1);
190 XSetForeground(disp, gc, BLACK);
191 XDrawCircle(disp, w, gc, xi, yi, size[ai]);
192 /* XSetForeground(disp,gc,WHITE);
193 XFillCircle(disp,w,gc,xi+4,yi-4,4); */
195 else if (bPlus)
197 XDrawLine(disp, w, gc, xi-MSIZE, yi, xi+MSIZE+1, yi);
198 XDrawLine(disp, w, gc, xi, yi-MSIZE, xi, yi+MSIZE+1);
200 else
202 XDrawLine(disp, w, gc, xi-1, yi, xi+1, yi);
207 /* Global variables */
208 static rvec gl_fbox, gl_hbox, gl_mhbox;
210 static void my_init_pbc(matrix box)
212 int i;
214 for (i = 0; (i < DIM); i++)
216 gl_fbox[i] = box[i][i];
217 gl_hbox[i] = gl_fbox[i]*0.5;
218 gl_mhbox[i] = -gl_hbox[i];
222 static bool local_pbc_dx(rvec x1, rvec x2)
224 int i;
225 real dx;
227 for (i = 0; (i < DIM); i++)
229 dx = x1[i]-x2[i];
230 if (dx > gl_hbox[i])
232 return false;
234 else if (dx <= gl_mhbox[i])
236 return false;
239 return true;
242 static void draw_bond(Display *disp, Window w, GC gc,
243 int ai, int aj, iv2 vec2[],
244 rvec x[], unsigned long col[], int size[], bool bBalls)
246 unsigned long ic, jc;
247 int xi, yi, xj, yj;
248 int xm, ym;
250 if (bBalls)
252 draw_atom(disp, w, gc, ai, vec2, col, size, true, false);
253 draw_atom(disp, w, gc, aj, vec2, col, size, true, false);
255 else
257 if (local_pbc_dx(x[ai], x[aj]))
259 ic = col[ai];
260 jc = col[aj];
261 xi = vec2[ai][XX];
262 yi = vec2[ai][YY];
263 xj = vec2[aj][XX];
264 yj = vec2[aj][YY];
266 if (ic != jc)
268 xm = (xi+xj) >> 1;
269 ym = (yi+yj) >> 1;
271 XSetForeground(disp, gc, ic);
272 XDrawLine(disp, w, gc, xi, yi, xm, ym);
273 XSetForeground(disp, gc, jc);
274 XDrawLine(disp, w, gc, xm, ym, xj, yj);
276 else
278 XSetForeground(disp, gc, ic);
279 XDrawLine(disp, w, gc, xi, yi, xj, yj);
285 int compare_obj(const void *a, const void *b)
287 t_object *oa, *ob;
288 real z;
290 oa = (t_object *)a;
291 ob = (t_object *)b;
293 z = oa->z-ob->z;
295 if (z < 0)
297 return 1;
299 else if (z > 0)
301 return -1;
303 else
305 return 0;
309 void z_fill(t_manager *man, real *zz)
311 t_object *obj;
312 int i;
314 for (i = 0, obj = man->obj; (i < man->nobj); i++, obj++)
316 switch (obj->eO)
318 case eOSingle:
319 obj->z = zz[obj->ai];
320 break;
321 case eOBond:
322 case eOHBond:
323 obj->z = (zz[obj->ai] + zz[obj->aj]) * 0.5;
324 break;
325 default:
326 break;
331 int filter_vis(t_manager *man)
333 int i, nobj, nvis, nhide;
334 int ai;
335 bool bAdd, *bVis;
336 t_object *obj;
337 t_object *newobj;
339 nobj = man->nobj;
340 snew(newobj, nobj);
341 obj = man->obj;
342 bVis = man->bVis;
343 nvis = 0;
344 nhide = nobj-1;
345 for (i = 0; (i < nobj); i++, obj++)
347 ai = obj->ai;
348 bAdd = bVis[ai];
349 if (bAdd)
351 if (obj->eO != eOSingle)
353 bAdd = bVis[obj->aj];
356 if (bAdd)
358 newobj[nvis++] = *obj;
360 else
362 newobj[nhide--] = *obj;
365 sfree(man->obj);
366 man->obj = newobj;
368 return nvis;
371 void draw_objects(Display *disp, Window w, GC gc, int nobj,
372 t_object objs[], iv2 vec2[], rvec x[],
373 unsigned long col[], int size[], bool bShowHydro, int bond_type,
374 bool bPlus)
376 bool bBalls;
377 int i;
378 t_object *obj;
380 bBalls = false;
381 switch (bond_type)
383 case eBThin:
384 XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound);
385 break;
386 case eBFat:
387 XSetLineAttributes(disp, gc, 3, LineSolid, CapNotLast, JoinRound);
388 break;
389 case eBVeryFat:
390 XSetLineAttributes(disp, gc, 5, LineSolid, CapNotLast, JoinRound);
391 break;
392 case eBSpheres:
393 bBalls = true;
394 bPlus = false;
395 break;
396 default:
397 gmx_fatal(FARGS, "Invalid bond_type selected: %d\n", bond_type);
398 break;
400 for (i = 0; (i < nobj); i++)
402 obj = &(objs[i]);
403 switch (obj->eO)
405 case eOSingle:
406 draw_atom(disp, w, gc, obj->ai, vec2, col, size, bBalls, bPlus);
407 break;
408 case eOBond:
409 draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
410 break;
411 case eOHBond:
412 if (bShowHydro)
414 draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
416 break;
417 default:
418 break;
421 XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound);
424 static void v4_to_iv2(vec4 x4, iv2 v2, int x0, int y0, real sx, real sy)
426 real inv_z;
428 inv_z = 1.0/x4[ZZ];
429 v2[XX] = x0+sx*x4[XX]*inv_z;
430 v2[YY] = y0-sy*x4[YY]*inv_z;
433 static void draw_box(t_x11 *x11, Window w, t_3dview *view, matrix box,
434 int x0, int y0, real sx, real sy, int boxtype)
436 rvec rect_tri[8] = {
437 { 0, 0, 0 }, { 1, 0, 0 }, { 1, 1, 0 }, { 0, 1, 0 },
438 { 0, 0, 1 }, { 1, 0, 1 }, { 1, 1, 1 }, { 0, 1, 1 }
440 int tr_bonds[12][2] = {
441 { 0, 1 }, { 1, 2 }, { 2, 3 }, { 3, 0 },
442 { 4, 5 }, { 5, 6 }, { 6, 7 }, { 7, 4 },
443 { 0, 4 }, { 1, 5 }, { 2, 6 }, { 3, 7 }
445 static int *edge = NULL;
446 int i, j, k, i0, i1;
447 rvec corner[NCUCEDGE], box_center;
448 vec4 x4;
449 iv2 vec2[NCUCEDGE];
451 calc_box_center(view->ecenter, box, box_center);
452 if (boxtype == esbTrunc)
454 calc_compact_unitcell_vertices(view->ecenter, box, corner);
455 if (edge == NULL)
457 edge = compact_unitcell_edges();
460 for (i = 0; (i < NCUCEDGE); i++)
462 gmx_mat4_transform_point(view->proj, corner[i], x4);
463 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
465 XSetForeground(x11->disp, x11->gc, YELLOW);
466 for (i = 0; i < NCUCEDGE; i++)
468 i0 = edge[2*i];
469 i1 = edge[2*i+1];
470 XDrawLine(x11->disp, w, x11->gc,
471 vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
474 else
476 if (boxtype == esbRect)
478 for (j = 0; (j < DIM); j++)
480 box_center[j] -= 0.5*box[j][j];
483 else
485 for (i = 0; (i < DIM); i++)
487 for (j = 0; (j < DIM); j++)
489 box_center[j] -= 0.5*box[i][j];
493 for (i = 0; (i < 8); i++)
495 clear_rvec(corner[i]);
496 for (j = 0; (j < DIM); j++)
498 if (boxtype == esbTri)
500 for (k = 0; (k < DIM); k++)
502 corner[i][k] += rect_tri[i][j]*box[j][k];
505 else
507 corner[i][j] = rect_tri[i][j]*box[j][j];
510 rvec_inc(corner[i], box_center);
511 gmx_mat4_transform_point(view->proj, corner[i], x4);
512 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
514 if (debug)
516 pr_rvecs(debug, 0, "box", box, DIM);
517 pr_rvecs(debug, 0, "corner", corner, 8);
519 XSetForeground(x11->disp, x11->gc, YELLOW);
520 for (i = 0; (i < 12); i++)
522 i0 = tr_bonds[i][0];
523 i1 = tr_bonds[i][1];
524 XDrawLine(x11->disp, w, x11->gc,
525 vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
530 void set_sizes(t_manager *man)
532 for (int i = 0; i < man->natom; i++)
534 if (man->bVis[i])
536 man->size[i] = 180*man->vdw[i];
541 void draw_mol(t_x11 *x11, t_manager *man)
543 static char tstr[2][20];
544 static int ntime = 0;
545 t_windata *win;
546 t_3dview *view;
547 t_molwin *mw;
548 int i, x0, y0, nvis;
549 iv2 *vec2;
550 real sx, sy;
551 vec4 x4;
553 if (!man->status)
555 return;
558 view = man->view;
559 mw = man->molw;
561 win = &(mw->wd);
563 vec2 = man->ix;
564 x0 = win->width/2;
565 y0 = win->height/2;
566 sx = win->width/2*view->sc_x;
567 sy = win->height/2*view->sc_y;
569 my_init_pbc(man->box);
571 for (i = 0; (i < man->natom); i++)
573 if (man->bVis[i])
575 gmx_mat4_transform_point(view->proj, man->x[i], x4);
576 man->zz[i] = x4[ZZ];
577 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
580 set_sizes(man);
582 z_fill (man, man->zz);
584 /* Start drawing */
585 XClearWindow(x11->disp, win->self);
587 /* Draw Time */
588 std::sprintf(tstr[ntime], "Time: %.3f ps", man->time);
589 if (std::strcmp(tstr[ntime], tstr[1-ntime]) != 0)
591 set_vbtime(x11, man->vbox, tstr[ntime]);
592 ntime = 1-ntime;
595 if (mw->boxtype != esbNone)
597 draw_box(x11, win->self, view, man->box, x0, y0, sx, sy, mw->boxtype);
600 /* Should sort on Z-Coordinates here! */
601 nvis = filter_vis(man);
602 if (nvis && man->bSort)
604 std::qsort(man->obj, nvis, sizeof(man->obj[0]), compare_obj);
607 /* Draw the objects */
608 draw_objects(x11->disp, win->self, x11->gc,
609 nvis, man->obj, man->ix, man->x, man->col, man->size,
610 mw->bShowHydrogen, mw->bond_type, man->bPlus);
612 /* Draw the labels */
613 XSetForeground(x11->disp, x11->gc, WHITE);
614 for (i = 0; (i < man->natom); i++)
616 if (man->bLabel[i] && man->bVis[i])
618 XDrawString(x11->disp, win->self, x11->gc, vec2[i][XX]+2, vec2[i][YY]-2,
619 man->szLab[i], std::strlen(man->szLab[i]));
623 XSetForeground(x11->disp, x11->gc, x11->fg);