21f502e29cb9c39a038fbb1931ec514c65ad3fec
[alexxy/gromacs.git] / src / programs / view / nmol.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-2013, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2019, 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 "nmol.h"
40
41 #include <cstdio>
42 #include <cstdlib>
43 #include <cstring>
44
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"
50
51 #include "3dview.h"
52 #include "buttons.h"
53 #include "manager.h"
54 #include "xutil.h"
55
56 #define MSIZE 4
57
58 static bool MWCallBack(t_x11* x11, XEvent* event, Window /*w*/, void* data)
59 {
60     t_molwin* mw;
61     Window    To;
62     XEvent    letter;
63
64     mw                          = static_cast<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)
72     {
73         case Expose:
74             /* Do Not draw anything, but signal parent instead, he will
75              * coordinate drawing.
76              */
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: break;
96     }
97     return false;
98 }
99
100 static void set_def(t_molwin* mw, int ePBC, matrix box)
101 {
102     mw->bShowHydrogen = true;
103     mw->bond_type     = eBFat;
104     mw->ePBC          = ePBC;
105     mw->boxtype       = esbRect;
106     mw->realbox       = TRICLINIC(box) ? esbTri : esbRect;
107 }
108
109 t_molwin* init_mw(t_x11*        x11,
110                   Window        Parent,
111                   int           x,
112                   int           y,
113                   int           width,
114                   int           height,
115                   unsigned long fg,
116                   unsigned long bg,
117                   int           ePBC,
118                   matrix        box)
119 {
120     t_molwin* mw;
121
122     snew(mw, 1);
123     set_def(mw, ePBC, box);
124
125     InitWin(&mw->wd, x, y, width, height, 1, "Mol Window");
126
127     mw->wd.Parent = Parent;
128     mw->wd.self   = XCreateSimpleWindow(x11->disp, Parent, x, y, width, height, 1, fg, bg);
129     x11->RegisterCallback(x11, mw->wd.self, Parent, MWCallBack, mw);
130     x11->SetInputMask(x11, mw->wd.self, ExposureMask | StructureNotifyMask | ButtonPressMask);
131     return mw;
132 }
133
134 void map_mw(t_x11* x11, t_molwin* mw)
135 {
136     XMapWindow(x11->disp, mw->wd.self);
137 }
138
139 bool toggle_hydrogen(t_x11* x11, t_molwin* mw)
140 {
141     mw->bShowHydrogen = !mw->bShowHydrogen;
142     ExposeWin(x11->disp, mw->wd.self);
143
144     return mw->bShowHydrogen;
145 }
146
147 void set_bond_type(t_x11* x11, t_molwin* mw, int bt)
148 {
149     if (bt != mw->bond_type)
150     {
151         mw->bond_type = bt;
152         ExposeWin(x11->disp, mw->wd.self);
153     }
154 }
155
156 void set_box_type(t_x11* x11, t_molwin* mw, int bt)
157 {
158 #ifdef DEBUG
159     std::fprintf(stderr, "mw->boxtype = %d, bt = %d\n", mw->boxtype, bt);
160 #endif
161     if (bt != mw->boxtype)
162     {
163         if ((bt == esbTrunc && mw->realbox == esbTri) || bt == esbTri || bt == esbNone)
164         {
165             mw->boxtype = bt;
166             ExposeWin(x11->disp, mw->wd.self);
167         }
168         else
169         {
170             std::fprintf(stderr, "Can not change rectangular box to truncated octahedron\n");
171         }
172     }
173 }
174
175 void done_mw(t_x11* x11, t_molwin* mw)
176 {
177     x11->UnRegisterCallback(x11, mw->wd.self);
178     sfree(mw);
179 }
180
181 static void
182 draw_atom(Display* disp, Window w, GC gc, int ai, iv2 vec2[], unsigned long col[], int size[], bool bBall, bool bPlus)
183 {
184     int xi, yi;
185
186     xi = vec2[ai][XX];
187     yi = vec2[ai][YY];
188     XSetForeground(disp, gc, col[ai]);
189     if (bBall)
190     {
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); */
196     }
197     else if (bPlus)
198     {
199         XDrawLine(disp, w, gc, xi - MSIZE, yi, xi + MSIZE + 1, yi);
200         XDrawLine(disp, w, gc, xi, yi - MSIZE, xi, yi + MSIZE + 1);
201     }
202     else
203     {
204         XDrawLine(disp, w, gc, xi - 1, yi, xi + 1, yi);
205     }
206 }
207
208 /* Global variables */
209 static rvec gl_fbox, gl_hbox, gl_mhbox;
210
211 static void my_init_pbc(matrix box)
212 {
213     int i;
214
215     for (i = 0; (i < DIM); i++)
216     {
217         gl_fbox[i]  = box[i][i];
218         gl_hbox[i]  = gl_fbox[i] * 0.5;
219         gl_mhbox[i] = -gl_hbox[i];
220     }
221 }
222
223 static bool local_pbc_dx(rvec x1, rvec x2)
224 {
225     int  i;
226     real dx;
227
228     for (i = 0; (i < DIM); i++)
229     {
230         dx = x1[i] - x2[i];
231         if (dx > gl_hbox[i])
232         {
233             return false;
234         }
235         else if (dx <= gl_mhbox[i])
236         {
237             return false;
238         }
239     }
240     return true;
241 }
242
243 static void
244 draw_bond(Display* disp, Window w, GC gc, int ai, int aj, iv2 vec2[], rvec x[], unsigned long col[], int size[], bool bBalls)
245 {
246     unsigned long ic, jc;
247     int           xi, yi, xj, yj;
248     int           xm, ym;
249
250     if (bBalls)
251     {
252         draw_atom(disp, w, gc, ai, vec2, col, size, true, false);
253         draw_atom(disp, w, gc, aj, vec2, col, size, true, false);
254     }
255     else
256     {
257         if (local_pbc_dx(x[ai], x[aj]))
258         {
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];
265
266             if (ic != jc)
267             {
268                 xm = (xi + xj) >> 1;
269                 ym = (yi + yj) >> 1;
270
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);
275             }
276             else
277             {
278                 XSetForeground(disp, gc, ic);
279                 XDrawLine(disp, w, gc, xi, yi, xj, yj);
280             }
281         }
282     }
283 }
284
285 int compare_obj(const void* a, const void* b)
286 {
287     const t_object *oa, *ob;
288     real            z;
289
290     oa = static_cast<const t_object*>(a);
291     ob = static_cast<const t_object*>(b);
292
293     z = oa->z - ob->z;
294
295     if (z < 0)
296     {
297         return 1;
298     }
299     else if (z > 0)
300     {
301         return -1;
302     }
303     else
304     {
305         return 0;
306     }
307 }
308
309 void z_fill(t_manager* man, real* zz)
310 {
311     t_object* obj;
312     int       i;
313
314     for (i = 0, obj = man->obj; (i < man->nobj); i++, obj++)
315     {
316         switch (obj->eO)
317         {
318             case eOSingle: obj->z = zz[obj->ai]; break;
319             case eOBond:
320             case eOHBond: obj->z = (zz[obj->ai] + zz[obj->aj]) * 0.5; break;
321             default: break;
322         }
323     }
324 }
325
326 int filter_vis(t_manager* man)
327 {
328     int       i, nobj, nvis, nhide;
329     int       ai;
330     bool      bAdd, *bVis;
331     t_object* obj;
332     t_object* newobj;
333
334     nobj = man->nobj;
335     snew(newobj, nobj);
336     obj   = man->obj;
337     bVis  = man->bVis;
338     nvis  = 0;
339     nhide = nobj - 1;
340     for (i = 0; (i < nobj); i++, obj++)
341     {
342         ai   = obj->ai;
343         bAdd = bVis[ai];
344         if (bAdd)
345         {
346             if (obj->eO != eOSingle)
347             {
348                 bAdd = bVis[obj->aj];
349             }
350         }
351         if (bAdd)
352         {
353             newobj[nvis++] = *obj;
354         }
355         else
356         {
357             newobj[nhide--] = *obj;
358         }
359     }
360     sfree(man->obj);
361     man->obj = newobj;
362
363     return nvis;
364 }
365
366 static void draw_objects(Display*      disp,
367                          Window        w,
368                          GC            gc,
369                          int           nobj,
370                          t_object      objs[],
371                          iv2           vec2[],
372                          rvec          x[],
373                          unsigned long col[],
374                          int           size[],
375                          bool          bShowHydro,
376                          int           bond_type,
377                          bool          bPlus)
378 {
379     bool      bBalls;
380     int       i;
381     t_object* obj;
382
383     bBalls = false;
384     switch (bond_type)
385     {
386         case eBThin: XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound); break;
387         case eBFat: XSetLineAttributes(disp, gc, 3, LineSolid, CapNotLast, JoinRound); break;
388         case eBVeryFat: XSetLineAttributes(disp, gc, 5, LineSolid, CapNotLast, JoinRound); break;
389         case eBSpheres:
390             bBalls = true;
391             bPlus  = false;
392             break;
393         default: gmx_fatal(FARGS, "Invalid bond_type selected: %d\n", bond_type); break;
394     }
395     for (i = 0; (i < nobj); i++)
396     {
397         obj = &(objs[i]);
398         switch (obj->eO)
399         {
400             case eOSingle: draw_atom(disp, w, gc, obj->ai, vec2, col, size, bBalls, bPlus); break;
401             case eOBond:
402                 draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
403                 break;
404             case eOHBond:
405                 if (bShowHydro)
406                 {
407                     draw_bond(disp, w, gc, obj->ai, obj->aj, vec2, x, col, size, bBalls);
408                 }
409                 break;
410             default: break;
411         }
412     }
413     XSetLineAttributes(disp, gc, 1, LineSolid, CapNotLast, JoinRound);
414 }
415
416 static void v4_to_iv2(vec4 x4, iv2 v2, int x0, int y0, real sx, real sy)
417 {
418     real inv_z;
419
420     inv_z  = 1.0 / x4[ZZ];
421     v2[XX] = x0 + sx * x4[XX] * inv_z;
422     v2[YY] = y0 - sy * x4[YY] * inv_z;
423 }
424
425 static void draw_box(t_x11* x11, Window w, t_3dview* view, matrix box, int x0, int y0, real sx, real sy, int boxtype)
426 {
427     rvec        rect_tri[8]     = { { 0, 0, 0 }, { 1, 0, 0 }, { 1, 1, 0 }, { 0, 1, 0 },
428                          { 0, 0, 1 }, { 1, 0, 1 }, { 1, 1, 1 }, { 0, 1, 1 } };
429     int         tr_bonds[12][2] = { { 0, 1 }, { 1, 2 }, { 2, 3 }, { 3, 0 }, { 4, 5 }, { 5, 6 },
430                             { 6, 7 }, { 7, 4 }, { 0, 4 }, { 1, 5 }, { 2, 6 }, { 3, 7 } };
431     static int* edge            = nullptr;
432     int         i, j, k, i0, i1;
433     rvec        corner[NCUCEDGE], box_center;
434     vec4        x4;
435     iv2         vec2[NCUCEDGE];
436
437     calc_box_center(view->ecenter, box, box_center);
438     if (boxtype == esbTrunc)
439     {
440         calc_compact_unitcell_vertices(view->ecenter, box, corner);
441         if (edge == nullptr)
442         {
443             edge = compact_unitcell_edges();
444         }
445
446         for (i = 0; (i < NCUCEDGE); i++)
447         {
448             gmx_mat4_transform_point(view->proj, corner[i], x4);
449             v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
450         }
451         XSetForeground(x11->disp, x11->gc, YELLOW);
452         for (i = 0; i < NCUCEDGE; i++)
453         {
454             i0 = edge[2 * i];
455             i1 = edge[2 * i + 1];
456             XDrawLine(x11->disp, w, x11->gc, vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
457         }
458     }
459     else
460     {
461         if (boxtype == esbRect)
462         {
463             for (j = 0; (j < DIM); j++)
464             {
465                 box_center[j] -= 0.5 * box[j][j];
466             }
467         }
468         else
469         {
470             for (i = 0; (i < DIM); i++)
471             {
472                 for (j = 0; (j < DIM); j++)
473                 {
474                     box_center[j] -= 0.5 * box[i][j];
475                 }
476             }
477         }
478         for (i = 0; (i < 8); i++)
479         {
480             clear_rvec(corner[i]);
481             for (j = 0; (j < DIM); j++)
482             {
483                 if (boxtype == esbTri)
484                 {
485                     for (k = 0; (k < DIM); k++)
486                     {
487                         corner[i][k] += rect_tri[i][j] * box[j][k];
488                     }
489                 }
490                 else
491                 {
492                     corner[i][j] = rect_tri[i][j] * box[j][j];
493                 }
494             }
495             rvec_inc(corner[i], box_center);
496             gmx_mat4_transform_point(view->proj, corner[i], x4);
497             v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
498         }
499         if (debug)
500         {
501             pr_rvecs(debug, 0, "box", box, DIM);
502             pr_rvecs(debug, 0, "corner", corner, 8);
503         }
504         XSetForeground(x11->disp, x11->gc, YELLOW);
505         for (i = 0; (i < 12); i++)
506         {
507             i0 = tr_bonds[i][0];
508             i1 = tr_bonds[i][1];
509             XDrawLine(x11->disp, w, x11->gc, vec2[i0][XX], vec2[i0][YY], vec2[i1][XX], vec2[i1][YY]);
510         }
511     }
512 }
513
514 void set_sizes(t_manager* man)
515 {
516     for (int i = 0; i < man->natom; i++)
517     {
518         if (man->bVis[i])
519         {
520             man->size[i] = 180 * man->vdw[i];
521         }
522     }
523 }
524
525 void draw_mol(t_x11* x11, t_manager* man)
526 {
527     static char tstr[2][20];
528     static int  ntime = 0;
529     t_windata*  win;
530     t_3dview*   view;
531     t_molwin*   mw;
532     int         i, x0, y0, nvis;
533     iv2*        vec2;
534     real        sx, sy;
535     vec4        x4;
536
537     if (!man->status)
538     {
539         return;
540     }
541
542     view = man->view;
543     mw   = man->molw;
544
545     win = &(mw->wd);
546
547     vec2 = man->ix;
548     x0   = win->width / 2;
549     y0   = win->height / 2;
550     sx   = win->width / 2 * view->sc_x;
551     sy   = win->height / 2 * view->sc_y;
552
553     my_init_pbc(man->box);
554
555     for (i = 0; (i < man->natom); i++)
556     {
557         if (man->bVis[i])
558         {
559             gmx_mat4_transform_point(view->proj, man->x[i], x4);
560             man->zz[i] = x4[ZZ];
561             v4_to_iv2(x4, vec2[i], x0, y0, sx, sy);
562         }
563     }
564     set_sizes(man);
565
566     z_fill(man, man->zz);
567
568     /* Start drawing */
569     XClearWindow(x11->disp, win->self);
570
571     /* Draw Time */
572     std::sprintf(tstr[ntime], "Time: %.3f ps", man->time);
573     if (std::strcmp(tstr[ntime], tstr[1 - ntime]) != 0)
574     {
575         set_vbtime(x11, man->vbox, tstr[ntime]);
576         ntime = 1 - ntime;
577     }
578
579     if (mw->boxtype != esbNone)
580     {
581         draw_box(x11, win->self, view, man->box, x0, y0, sx, sy, mw->boxtype);
582     }
583
584     /* Should sort on Z-Coordinates here! */
585     nvis = filter_vis(man);
586     if (nvis && man->bSort)
587     {
588         std::qsort(man->obj, nvis, sizeof(man->obj[0]), compare_obj);
589     }
590
591     /* Draw the objects */
592     draw_objects(x11->disp, win->self, x11->gc, nvis, man->obj, man->ix, man->x, man->col,
593                  man->size, mw->bShowHydrogen, mw->bond_type, man->bPlus);
594
595     /* Draw the labels */
596     XSetForeground(x11->disp, x11->gc, WHITE);
597     for (i = 0; (i < man->natom); i++)
598     {
599         if (man->bLabel[i] && man->bVis[i])
600         {
601             XDrawString(x11->disp, win->self, x11->gc, vec2[i][XX] + 2, vec2[i][YY] - 2,
602                         man->szLab[i], std::strlen(man->szLab[i]));
603         }
604     }
605
606     XSetForeground(x11->disp, x11->gc, x11->fg);
607 }