SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_bundle.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) 2013,2014,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42 #include <vector>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/rmpbc.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/trajectory/trajectoryframe.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61
62 #define MAX_ENDS 3
63
64 typedef struct
65 {
66     int   n;
67     int   nend;
68     rvec* end[MAX_ENDS];
69     rvec* mid;
70     rvec* dir;
71     real* len;
72 } t_bundle;
73
74 static void rotate_ends(t_bundle* bun, rvec axis, int c0, int c1)
75 {
76     int  end, i;
77     rvec ax, tmp;
78
79     unitv(axis, ax);
80     for (end = 0; end < bun->nend; end++)
81     {
82         for (i = 0; i < bun->n; i++)
83         {
84             copy_rvec(bun->end[end][i], tmp);
85             bun->end[end][i][c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
86             bun->end[end][i][c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
87         }
88     }
89     copy_rvec(axis, tmp);
90     axis[c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
91     axis[c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
92 }
93
94 static void calc_axes(rvec x[], t_atom atom[], const int gnx[], int* index[], gmx_bool bRot, t_bundle* bun)
95 {
96     int   end, i, div, d;
97     real *mtot, m;
98     rvec  axis[MAX_ENDS], cent;
99
100     snew(mtot, bun->n);
101
102     clear_rvec(axis[0]);
103     clear_rvec(axis[1]);
104
105     for (end = 0; end < bun->nend; end++)
106     {
107         for (i = 0; i < bun->n; i++)
108         {
109             clear_rvec(bun->end[end][i]);
110             mtot[i] = 0;
111         }
112         div = gnx[end] / bun->n;
113         for (i = 0; i < gnx[end]; i++)
114         {
115             m = atom[index[end][i]].m;
116             for (d = 0; d < DIM; d++)
117             {
118                 bun->end[end][i / div][d] += m * x[index[end][i]][d];
119             }
120             mtot[i / div] += m;
121         }
122         clear_rvec(axis[end]);
123         for (i = 0; i < bun->n; i++)
124         {
125             svmul(1.0 / mtot[i], bun->end[end][i], bun->end[end][i]);
126             rvec_inc(axis[end], bun->end[end][i]);
127         }
128         svmul(1.0 / bun->n, axis[end], axis[end]);
129     }
130     sfree(mtot);
131
132     rvec_add(axis[0], axis[1], cent);
133     svmul(0.5, cent, cent);
134     /* center the bundle on the origin */
135     for (end = 0; end < bun->nend; end++)
136     {
137         rvec_dec(axis[end], cent);
138         for (i = 0; i < bun->n; i++)
139         {
140             rvec_dec(bun->end[end][i], cent);
141         }
142     }
143     if (bRot)
144     {
145         /* rotate the axis parallel to the z-axis */
146         rotate_ends(bun, axis[0], YY, ZZ);
147         rotate_ends(bun, axis[0], XX, ZZ);
148     }
149     for (i = 0; i < bun->n; i++)
150     {
151         rvec_add(bun->end[0][i], bun->end[1][i], bun->mid[i]);
152         svmul(0.5, bun->mid[i], bun->mid[i]);
153         rvec_sub(bun->end[0][i], bun->end[1][i], bun->dir[i]);
154         bun->len[i] = norm(bun->dir[i]);
155         unitv(bun->dir[i], bun->dir[i]);
156     }
157 }
158
159 static void dump_axes(t_trxstatus* status, t_trxframe* fr, t_atoms* outat, t_bundle* bun)
160 {
161     t_trxframe                    frout;
162     static std::vector<gmx::RVec> xout;
163     int                           i;
164
165     GMX_ASSERT(outat->nr >= bun->n, "");
166     if (xout.empty())
167     {
168         xout.resize(outat->nr);
169     }
170
171     for (i = 0; i < bun->n; i++)
172     {
173         copy_rvec(bun->end[0][i], xout[3 * i]);
174         if (bun->nend >= 3)
175         {
176             copy_rvec(bun->end[2][i], xout[3 * i + 1]);
177         }
178         else
179         {
180             copy_rvec(bun->mid[i], xout[3 * i + 1]);
181         }
182         copy_rvec(bun->end[1][i], xout[3 * i + 2]);
183     }
184     frout        = *fr;
185     frout.bV     = FALSE;
186     frout.bF     = FALSE;
187     frout.bBox   = FALSE;
188     frout.bAtoms = TRUE;
189     frout.natoms = outat->nr;
190     frout.atoms  = outat;
191     frout.x      = as_rvec_array(xout.data());
192     write_trxframe(status, &frout, nullptr);
193 }
194
195 int gmx_bundle(int argc, char* argv[])
196 {
197     const char* desc[] = {
198         "[THISMODULE] analyzes bundles of axes. The axes can be for instance",
199         "helix axes. The program reads two index groups and divides both",
200         "of them in [TT]-na[tt] parts. The centers of mass of these parts",
201         "define the tops and bottoms of the axes.",
202         "Several quantities are written to file:",
203         "the axis length, the distance and the z-shift of the axis mid-points",
204         "with respect to the average center of all axes, the total tilt,",
205         "the radial tilt and the lateral tilt with respect to the average axis.",
206         "[PAR]",
207         "With options [TT]-ok[tt], [TT]-okr[tt] and [TT]-okl[tt] the total,",
208         "radial and lateral kinks of the axes are plotted. An extra index",
209         "group of kink atoms is required, which is also divided into [TT]-na[tt]",
210         "parts. The kink angle is defined as the angle between the kink-top and",
211         "the bottom-kink vectors.",
212         "[PAR]",
213         "With option [TT]-oa[tt] the top, mid (or kink when [TT]-ok[tt] is set)",
214         "and bottom points of each axis",
215         "are written to a [REF].pdb[ref] file each frame. The residue numbers correspond",
216         "to the axis numbers. When viewing this file with Rasmol, use the",
217         "command line option [TT]-nmrpdb[tt], and type [TT]set axis true[tt] to",
218         "display the reference axis."
219     };
220     static int      n    = 0;
221     static gmx_bool bZ   = FALSE;
222     t_pargs         pa[] = { { "-na", FALSE, etINT, { &n }, "Number of axes" },
223                      { "-z",
224                        FALSE,
225                        etBOOL,
226                        { &bZ },
227                        "Use the [IT]z[it]-axis as reference instead of the average axis" } };
228     FILE *          flen, *fdist, *fz, *ftilt, *ftiltr, *ftiltl;
229     FILE *          fkink = nullptr, *fkinkr = nullptr, *fkinkl = nullptr;
230     t_trxstatus*    status;
231     t_trxstatus*    fpdb;
232     t_topology      top;
233     PbcType         pbcType;
234     rvec*           xtop;
235     matrix          box;
236     t_trxframe      fr;
237     t_atoms         outatoms;
238     real            t, comp;
239     char*           grpname[MAX_ENDS];
240     /* FIXME: The constness should not be cast away */
241     char *            anm = const_cast<char*>("CA"), *rnm = const_cast<char*>("GLY");
242     int               i, gnx[MAX_ENDS];
243     int*              index[MAX_ENDS];
244     t_bundle          bun;
245     gmx_bool          bKink;
246     rvec              va, vb, vc, vr, vl;
247     gmx_output_env_t* oenv;
248     gmx_rmpbc_t       gpbc = nullptr;
249
250     t_filenm fnm[] = {
251         { efTRX, "-f", nullptr, ffREAD },        { efTPS, nullptr, nullptr, ffREAD },
252         { efNDX, nullptr, nullptr, ffOPTRD },    { efXVG, "-ol", "bun_len", ffWRITE },
253         { efXVG, "-od", "bun_dist", ffWRITE },   { efXVG, "-oz", "bun_z", ffWRITE },
254         { efXVG, "-ot", "bun_tilt", ffWRITE },   { efXVG, "-otr", "bun_tiltr", ffWRITE },
255         { efXVG, "-otl", "bun_tiltl", ffWRITE }, { efXVG, "-ok", "bun_kink", ffOPTWR },
256         { efXVG, "-okr", "bun_kinkr", ffOPTWR }, { efXVG, "-okl", "bun_kinkl", ffOPTWR },
257         { efPDB, "-oa", "axes", ffOPTWR }
258     };
259 #define NFILE asize(fnm)
260
261     if (!parse_common_args(
262                 &argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
263     {
264         return 0;
265     }
266
267     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &xtop, nullptr, box, TRUE);
268
269     bKink = opt2bSet("-ok", NFILE, fnm) || opt2bSet("-okr", NFILE, fnm) || opt2bSet("-okl", NFILE, fnm);
270     if (bKink)
271     {
272         bun.nend = 3;
273     }
274     else
275     {
276         bun.nend = 2;
277     }
278
279     fprintf(stderr, "Select a group of top and a group of bottom ");
280     if (bKink)
281     {
282         fprintf(stderr, "and a group of kink ");
283     }
284     fprintf(stderr, "atoms\n");
285     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bun.nend, gnx, index, grpname);
286
287     if (n <= 0 || gnx[0] % n || gnx[1] % n || (bKink && gnx[2] % n))
288     {
289         gmx_fatal(FARGS, "The size of one of your index groups is not a multiple of n");
290     }
291     bun.n = n;
292     snew(bun.end[0], n);
293     snew(bun.end[1], n);
294     if (bKink)
295     {
296         snew(bun.end[2], n);
297     }
298     snew(bun.mid, n);
299     snew(bun.dir, n);
300     snew(bun.len, n);
301
302     flen = xvgropen(
303             opt2fn("-ol", NFILE, fnm), "Axis lengths", output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
304     fdist = xvgropen(opt2fn("-od", NFILE, fnm),
305                      "Distance of axis centers",
306                      output_env_get_xvgr_tlabel(oenv),
307                      "(nm)",
308                      oenv);
309     fz    = xvgropen(opt2fn("-oz", NFILE, fnm),
310                   "Z-shift of axis centers",
311                   output_env_get_xvgr_tlabel(oenv),
312                   "(nm)",
313                   oenv);
314     ftilt = xvgropen(
315             opt2fn("-ot", NFILE, fnm), "Axis tilts", output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
316     ftiltr = xvgropen(opt2fn("-otr", NFILE, fnm),
317                       "Radial axis tilts",
318                       output_env_get_xvgr_tlabel(oenv),
319                       "(degrees)",
320                       oenv);
321     ftiltl = xvgropen(opt2fn("-otl", NFILE, fnm),
322                       "Lateral axis tilts",
323                       output_env_get_xvgr_tlabel(oenv),
324                       "(degrees)",
325                       oenv);
326
327     if (bKink)
328     {
329         fkink  = xvgropen(opt2fn("-ok", NFILE, fnm),
330                          "Kink angles",
331                          output_env_get_xvgr_tlabel(oenv),
332                          "(degrees)",
333                          oenv);
334         fkinkr = xvgropen(opt2fn("-okr", NFILE, fnm),
335                           "Radial kink angles",
336                           output_env_get_xvgr_tlabel(oenv),
337                           "(degrees)",
338                           oenv);
339         if (output_env_get_print_xvgr_codes(oenv))
340         {
341             fprintf(fkinkr, "@ subtitle \"+ = ) (   - = ( )\"\n");
342         }
343         fkinkl = xvgropen(opt2fn("-okl", NFILE, fnm),
344                           "Lateral kink angles",
345                           output_env_get_xvgr_tlabel(oenv),
346                           "(degrees)",
347                           oenv);
348     }
349
350     if (opt2bSet("-oa", NFILE, fnm))
351     {
352         init_t_atoms(&outatoms, 3 * n, FALSE);
353         outatoms.nr = 3 * n;
354         for (i = 0; i < 3 * n; i++)
355         {
356             outatoms.atomname[i]         = &anm;
357             outatoms.atom[i].resind      = i / 3;
358             outatoms.resinfo[i / 3].name = &rnm;
359             outatoms.resinfo[i / 3].nr   = i / 3 + 1;
360             outatoms.resinfo[i / 3].ic   = ' ';
361         }
362         fpdb = open_trx(opt2fn("-oa", NFILE, fnm), "w");
363     }
364     else
365     {
366         fpdb = nullptr;
367     }
368
369     read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, TRX_NEED_X);
370     gpbc = gmx_rmpbc_init(&top.idef, pbcType, fr.natoms);
371
372     do
373     {
374         gmx_rmpbc_trxfr(gpbc, &fr);
375         calc_axes(fr.x, top.atoms.atom, gnx, index, !bZ, &bun);
376         t = output_env_conv_time(oenv, fr.time);
377         fprintf(flen, " %10g", t);
378         fprintf(fdist, " %10g", t);
379         fprintf(fz, " %10g", t);
380         fprintf(ftilt, " %10g", t);
381         fprintf(ftiltr, " %10g", t);
382         fprintf(ftiltl, " %10g", t);
383         if (bKink)
384         {
385             fprintf(fkink, " %10g", t);
386             fprintf(fkinkr, " %10g", t);
387             fprintf(fkinkl, " %10g", t);
388         }
389
390         for (i = 0; i < bun.n; i++)
391         {
392             fprintf(flen, " %6g", bun.len[i]);
393             fprintf(fdist, " %6g", norm(bun.mid[i]));
394             fprintf(fz, " %6g", bun.mid[i][ZZ]);
395             fprintf(ftilt, " %6g", gmx::c_rad2Deg * acos(bun.dir[i][ZZ]));
396             comp = bun.mid[i][XX] * bun.dir[i][XX] + bun.mid[i][YY] * bun.dir[i][YY];
397             fprintf(ftiltr, " %6g", gmx::c_rad2Deg * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
398             comp = bun.mid[i][YY] * bun.dir[i][XX] - bun.mid[i][XX] * bun.dir[i][YY];
399             fprintf(ftiltl, " %6g", gmx::c_rad2Deg * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
400             if (bKink)
401             {
402                 rvec_sub(bun.end[0][i], bun.end[2][i], va);
403                 rvec_sub(bun.end[2][i], bun.end[1][i], vb);
404                 unitv(va, va);
405                 unitv(vb, vb);
406                 fprintf(fkink, " %6g", gmx::c_rad2Deg * acos(iprod(va, vb)));
407                 cprod(va, vb, vc);
408                 copy_rvec(bun.mid[i], vr);
409                 vr[ZZ] = 0;
410                 unitv(vr, vr);
411                 fprintf(fkinkr, " %6g", gmx::c_rad2Deg * std::asin(iprod(vc, vr)));
412                 vl[XX] = vr[YY];
413                 vl[YY] = -vr[XX];
414                 vl[ZZ] = 0;
415                 fprintf(fkinkl, " %6g", gmx::c_rad2Deg * std::asin(iprod(vc, vl)));
416             }
417         }
418         fprintf(flen, "\n");
419         fprintf(fdist, "\n");
420         fprintf(fz, "\n");
421         fprintf(ftilt, "\n");
422         fprintf(ftiltr, "\n");
423         fprintf(ftiltl, "\n");
424         if (bKink)
425         {
426             fprintf(fkink, "\n");
427             fprintf(fkinkr, "\n");
428             fprintf(fkinkl, "\n");
429         }
430         if (fpdb)
431         {
432             dump_axes(fpdb, &fr, &outatoms, &bun);
433         }
434     } while (read_next_frame(oenv, status, &fr));
435     gmx_rmpbc_done(gpbc);
436
437     close_trx(status);
438
439     if (fpdb)
440     {
441         close_trx(fpdb);
442     }
443     xvgrclose(flen);
444     xvgrclose(fdist);
445     xvgrclose(fz);
446     xvgrclose(ftilt);
447     xvgrclose(ftiltr);
448     xvgrclose(ftiltl);
449     if (bKink)
450     {
451         xvgrclose(fkink);
452         xvgrclose(fkinkr);
453         xvgrclose(fkinkl);
454     }
455
456     return 0;
457 }