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, 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.
50 #include "gromacs/math/vec.h"
54 #include "gromacs/math/3dview.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/smalloc.h"
60 static bool MWCallBack(t_x11 *x11, XEvent *event, Window /*w*/, void *data)
66 mw = (t_molwin *)data;
68 letter.type = ClientMessage;
69 letter.xclient.display = x11->disp;
70 letter.xclient.window = To;
71 letter.xclient.message_type = 0;
72 letter.xclient.format = 32;
76 /* Do Not draw anything, but signal parent instead, he will
79 letter.xclient.data.l[0] = IDDRAWMOL;
80 letter.xclient.data.l[1] = Button1;
81 XSendEvent(x11->disp, To, True, 0, &letter);
85 printf("Molwindow: Buttonpress\n");
87 letter.xclient.data.l[0] = IDLABEL;
88 letter.xclient.data.l[1] = (long)event->xbutton.button;
89 letter.xclient.data.l[2] = event->xbutton.x;
90 letter.xclient.data.l[3] = event->xbutton.y;
91 XSendEvent(x11->disp, To, True, 0, &letter);
94 mw->wd.width = event->xconfigure.width;
95 mw->wd.height = event->xconfigure.height;
103 void set_def (t_molwin *mw, int ePBC, matrix box)
105 mw->bShowHydrogen = true;
106 mw->bond_type = eBFat;
108 mw->boxtype = esbRect;
109 mw->realbox = TRICLINIC(box) ? esbTri : esbRect;
112 t_molwin *init_mw(t_x11 *x11, Window Parent,
113 int x, int y, int width, int height,
114 unsigned long fg, unsigned long bg,
115 int ePBC, matrix box)
120 set_def(mw, ePBC, box);
122 InitWin(&mw->wd, x, y, width, height, 1, "Mol Window");
124 mw->wd.Parent = Parent;
125 mw->wd.self = XCreateSimpleWindow(x11->disp, Parent, x, y, width, height, 1, fg, bg);
126 x11->RegisterCallback(x11, mw->wd.self, Parent, MWCallBack, mw);
127 x11->SetInputMask(x11, mw->wd.self,
128 ExposureMask | StructureNotifyMask |
133 void map_mw(t_x11 *x11, t_molwin *mw)
135 XMapWindow(x11->disp, mw->wd.self);
138 bool toggle_hydrogen(t_x11 *x11, t_molwin *mw)
140 mw->bShowHydrogen = !mw->bShowHydrogen;
141 ExposeWin(x11->disp, mw->wd.self);
143 return mw->bShowHydrogen;
146 void set_bond_type(t_x11 *x11, t_molwin *mw, int bt)
148 if (bt != mw->bond_type)
151 ExposeWin(x11->disp, mw->wd.self);
155 void set_box_type (t_x11 *x11, t_molwin *mw, int bt)
158 fprintf(stderr, "mw->boxtype = %d, bt = %d\n", mw->boxtype, bt);
160 if (bt != mw->boxtype)
162 if ((bt == esbTrunc && mw->realbox == esbTri) || bt == esbTri || bt == esbNone)
165 ExposeWin(x11->disp, mw->wd.self);
169 fprintf(stderr, "Can not change rectangular box to truncated octahedron\n");
174 void done_mw(t_x11 *x11, t_molwin *mw)
176 x11->UnRegisterCallback(x11, mw->wd.self);
180 static void draw_atom(Display *disp, Window w, GC gc,
181 atom_id ai, iv2 vec2[], unsigned long col[], int size[],
182 bool bBall, bool bPlus)
188 XSetForeground(disp, gc, col[ai]);
191 XFillCircle(disp, w, gc, xi, yi, size[ai]-1);
192 XSetForeground(disp, gc, BLACK);
193 XDrawCircle(disp, w, gc, xi, yi, size[ai]);
194 /* XSetForeground(disp,gc,WHITE);
195 XFillCircle(disp,w,gc,xi+4,yi-4,4); */
199 XDrawLine(disp, w, gc, xi-MSIZE, yi, xi+MSIZE+1, yi);
200 XDrawLine(disp, w, gc, xi, yi-MSIZE, xi, yi+MSIZE+1);
204 XDrawLine(disp, w, gc, xi-1, yi, xi+1, yi);
209 /* Global variables */
210 static rvec gl_fbox, gl_hbox, gl_mhbox;
212 static void my_init_pbc(matrix box)
216 for (i = 0; (i < DIM); i++)
218 gl_fbox[i] = box[i][i];
219 gl_hbox[i] = gl_fbox[i]*0.5;
220 gl_mhbox[i] = -gl_hbox[i];
224 static bool local_pbc_dx(rvec x1, rvec x2)
229 for (i = 0; (i < DIM); i++)
236 else if (dx <= gl_mhbox[i])
244 static void draw_bond(Display *disp, Window w, GC gc,
245 atom_id ai, atom_id aj, iv2 vec2[],
246 rvec x[], unsigned long col[], int size[], bool bBalls)
248 unsigned long ic, jc;
254 draw_atom(disp, w, gc, ai, vec2, col, size, true, false);
255 draw_atom(disp, w, gc, aj, vec2, col, size, true, false);
259 if (local_pbc_dx(x[ai], x[aj]))
273 XSetForeground(disp, gc, ic);
274 XDrawLine(disp, w, gc, xi, yi, xm, ym);
275 XSetForeground(disp, gc, jc);
276 XDrawLine(disp, w, gc, xm, ym, xj, yj);
280 XSetForeground(disp, gc, ic);
281 XDrawLine(disp, w, gc, xi, yi, xj, yj);
287 int compare_obj(const void *a, const void *b)
311 void z_fill(t_manager *man, real *zz)
316 for (i = 0, obj = man->obj; (i < man->nobj); i++, obj++)
321 obj->z = zz[obj->ai];
325 obj->z = (zz[obj->ai] + zz[obj->aj]) * 0.5;
333 int filter_vis(t_manager *man)
335 int i, nobj, nvis, nhide;
347 for (i = 0; (i < nobj); i++, obj++)
353 if (obj->eO != eOSingle)
355 bAdd = bVis[obj->aj];
360 newobj[nvis++] = *obj;
364 newobj[nhide--] = *obj;
373 void draw_objects(Display *disp, Window w, GC gc, int nobj,
374 t_object objs[], iv2 vec2[], rvec x[],
375 unsigned long col[], int size[], bool bShowHydro, int bond_type,
386 XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound);
389 XSetLineAttributes(disp, gc, 3, LineSolid, CapNotLast, JoinRound);
392 XSetLineAttributes(disp, gc, 5, LineSolid, CapNotLast, JoinRound);
399 gmx_fatal(FARGS, "Invalid bond_type selected: %d\n", bond_type);
402 for (i = 0; (i < nobj); i++)
408 draw_atom(disp, w, gc, obj->ai, vec2, col, size, bBalls, bPlus);
411 draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
416 draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
423 XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound);
426 static void v4_to_iv2(vec4 x4, iv2 v2, int x0, int y0, real sx, real sy)
431 v2[XX] = x0+sx*x4[XX]*inv_z;
432 v2[YY] = y0-sy*x4[YY]*inv_z;
435 static void draw_box(t_x11 *x11, Window w, t_3dview *view, matrix box,
436 int x0, int y0, real sx, real sy, int boxtype)
439 { 0, 0, 0 }, { 1, 0, 0 }, { 1, 1, 0 }, { 0, 1, 0 },
440 { 0, 0, 1 }, { 1, 0, 1 }, { 1, 1, 1 }, { 0, 1, 1 }
442 int tr_bonds[12][2] = {
443 { 0, 1 }, { 1, 2 }, { 2, 3 }, { 3, 0 },
444 { 4, 5 }, { 5, 6 }, { 6, 7 }, { 7, 4 },
445 { 0, 4 }, { 1, 5 }, { 2, 6 }, { 3, 7 }
447 static int *edge = NULL;
449 rvec corner[NCUCEDGE], box_center;
453 calc_box_center(view->ecenter, box, box_center);
454 if (boxtype == esbTrunc)
456 calc_compact_unitcell_vertices(view->ecenter, box, corner);
459 edge = compact_unitcell_edges();
462 for (i = 0; (i < NCUCEDGE); i++)
464 gmx_mat4_transform_point(view->proj, corner[i], x4);
465 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
467 XSetForeground(x11->disp, x11->gc, YELLOW);
468 for (i = 0; i < NCUCEDGE; i++)
472 XDrawLine(x11->disp, w, x11->gc,
473 vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
478 if (boxtype == esbRect)
480 for (j = 0; (j < DIM); j++)
482 box_center[j] -= 0.5*box[j][j];
487 for (i = 0; (i < DIM); i++)
489 for (j = 0; (j < DIM); j++)
491 box_center[j] -= 0.5*box[i][j];
495 for (i = 0; (i < 8); i++)
497 clear_rvec(corner[i]);
498 for (j = 0; (j < DIM); j++)
500 if (boxtype == esbTri)
502 for (k = 0; (k < DIM); k++)
504 corner[i][k] += rect_tri[i][j]*box[j][k];
509 corner[i][j] = rect_tri[i][j]*box[j][j];
512 rvec_inc(corner[i], box_center);
513 gmx_mat4_transform_point(view->proj, corner[i], x4);
514 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
518 pr_rvecs(debug, 0, "box", box, DIM);
519 pr_rvecs(debug, 0, "corner", corner, 8);
521 XSetForeground(x11->disp, x11->gc, YELLOW);
522 for (i = 0; (i < 12); i++)
526 XDrawLine(x11->disp, w, x11->gc,
527 vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
532 void set_sizes(t_manager *man)
534 for (int i = 0; i < man->natom; i++)
538 man->size[i] = 180*man->vdw[i];
543 void draw_mol(t_x11 *x11, t_manager *man)
545 static char tstr[2][20];
546 static int ntime = 0;
568 sx = win->width/2*view->sc_x;
569 sy = win->height/2*view->sc_y;
571 my_init_pbc(man->box);
573 for (i = 0; (i < man->natom); i++)
577 gmx_mat4_transform_point(view->proj, man->x[i], x4);
579 v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
584 z_fill (man, man->zz);
587 XClearWindow(x11->disp, win->self);
590 sprintf(tstr[ntime], "Time: %.3f ps", man->time);
591 if (strcmp(tstr[ntime], tstr[1-ntime]) != 0)
593 set_vbtime(x11, man->vbox, tstr[ntime]);
597 if (mw->boxtype != esbNone)
599 draw_box(x11, win->self, view, man->box, x0, y0, sx, sy, mw->boxtype);
602 /* Should sort on Z-Coordinates here! */
603 nvis = filter_vis(man);
604 if (nvis && man->bSort)
606 qsort(man->obj, nvis, sizeof(man->obj[0]), compare_obj);
609 /* Draw the objects */
610 draw_objects(x11->disp, win->self, x11->gc,
611 nvis, man->obj, man->ix, man->x, man->col, man->size,
612 mw->bShowHydrogen, mw->bond_type, man->bPlus);
614 /* Draw the labels */
615 XSetForeground(x11->disp, x11->gc, WHITE);
616 for (i = 0; (i < man->natom); i++)
618 if (man->bLabel[i] && man->bVis[i])
620 XDrawString(x11->disp, win->self, x11->gc, vec2[i][XX]+2, vec2[i][YY]-2,
621 man->szLab[i], strlen(man->szLab[i]));
625 XSetForeground(x11->disp, x11->gc, x11->fg);