Move mtop_util.* and topsort.* to topology/
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_nmeig.c
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <string.h>
43
44 #include "typedefs.h"
45 #include "macros.h"
46 #include "copyrite.h"
47 #include "gromacs/utility/futil.h"
48 #include "index.h"
49 #include "gstat.h"
50 #include "txtdump.h"
51 #include "eigio.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/math/units.h"
54 #include "gmx_ana.h"
55
56 #include "gromacs/commandline/pargs.h"
57 #include "gromacs/fileio/mtxio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/linearalgebra/eigensolver.h"
60 #include "gromacs/linearalgebra/sparsematrix.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/utility/smalloc.h"
63
64 static double cv_corr(double nu, double T)
65 {
66     double x  = PLANCK*nu/(BOLTZ*T);
67     double ex = exp(x);
68
69     if (nu <= 0)
70     {
71         return BOLTZ*KILO;
72     }
73     else
74     {
75         return BOLTZ*KILO*(ex*sqr(x)/sqr(ex-1) - 1);
76     }
77 }
78
79 static double u_corr(double nu, double T)
80 {
81     double x  = PLANCK*nu/(BOLTZ*T);
82     double ex = exp(x);
83
84     if (nu <= 0)
85     {
86         return BOLTZ*T;
87     }
88     else
89     {
90         return BOLTZ*T*(0.5*x - 1 + x/(ex-1));
91     }
92 }
93
94 static int get_nharm_mt(gmx_moltype_t *mt)
95 {
96     static int harm_func[] = { F_BONDS };
97     int        i, ft, nh;
98
99     nh = 0;
100     for (i = 0; (i < asize(harm_func)); i++)
101     {
102         ft  = harm_func[i];
103         nh += mt->ilist[ft].nr/(interaction_function[ft].nratoms+1);
104     }
105     return nh;
106 }
107
108 static int get_nvsite_mt(gmx_moltype_t *mt)
109 {
110     static int vs_func[] = {
111         F_VSITE2, F_VSITE3, F_VSITE3FD, F_VSITE3FAD,
112         F_VSITE3OUT, F_VSITE4FD, F_VSITE4FDN, F_VSITEN
113     };
114     int        i, ft, nh;
115
116     nh = 0;
117     for (i = 0; (i < asize(vs_func)); i++)
118     {
119         ft  = vs_func[i];
120         nh += mt->ilist[ft].nr/(interaction_function[ft].nratoms+1);
121     }
122     return nh;
123 }
124
125 static int get_nharm(gmx_mtop_t *mtop, int *nvsites)
126 {
127     int j, mt, nh, nv;
128
129     nh = 0;
130     nv = 0;
131     for (j = 0; (j < mtop->nmolblock); j++)
132     {
133         mt  = mtop->molblock[j].type;
134         nh += mtop->molblock[j].nmol * get_nharm_mt(&(mtop->moltype[mt]));
135         nv += mtop->molblock[j].nmol * get_nvsite_mt(&(mtop->moltype[mt]));
136     }
137     *nvsites = nv;
138     return nh;
139 }
140
141 static void
142 nma_full_hessian(real     *           hess,
143                  int                  ndim,
144                  gmx_bool             bM,
145                  t_topology     *     top,
146                  int                  begin,
147                  int                  end,
148                  real     *           eigenvalues,
149                  real     *           eigenvectors)
150 {
151     int  i, j, k, l;
152     real mass_fac, rdum;
153     int  natoms;
154
155     natoms = top->atoms.nr;
156
157     /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
158
159     if (bM)
160     {
161         for (i = 0; (i < natoms); i++)
162         {
163             for (j = 0; (j < DIM); j++)
164             {
165                 for (k = 0; (k < natoms); k++)
166                 {
167                     mass_fac = gmx_invsqrt(top->atoms.atom[i].m*top->atoms.atom[k].m);
168                     for (l = 0; (l < DIM); l++)
169                     {
170                         hess[(i*DIM+j)*ndim+k*DIM+l] *= mass_fac;
171                     }
172                 }
173             }
174         }
175     }
176
177     /* call diagonalization routine. */
178
179     fprintf(stderr, "\nDiagonalizing to find vectors %d through %d...\n", begin, end);
180     fflush(stderr);
181
182     eigensolver(hess, ndim, begin-1, end-1, eigenvalues, eigenvectors);
183
184     /* And scale the output eigenvectors */
185     if (bM && eigenvectors != NULL)
186     {
187         for (i = 0; i < (end-begin+1); i++)
188         {
189             for (j = 0; j < natoms; j++)
190             {
191                 mass_fac = gmx_invsqrt(top->atoms.atom[j].m);
192                 for (k = 0; (k < DIM); k++)
193                 {
194                     eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
195                 }
196             }
197         }
198     }
199 }
200
201
202
203 static void
204 nma_sparse_hessian(gmx_sparsematrix_t     *     sparse_hessian,
205                    gmx_bool                     bM,
206                    t_topology     *             top,
207                    int                          neig,
208                    real     *                   eigenvalues,
209                    real     *                   eigenvectors)
210 {
211     int  i, j, k;
212     int  row, col;
213     real mass_fac;
214     int  iatom, katom;
215     int  natoms;
216     int  ndim;
217
218     natoms = top->atoms.nr;
219     ndim   = DIM*natoms;
220
221     /* Cannot check symmetry since we only store half matrix */
222     /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
223
224     if (bM)
225     {
226         for (iatom = 0; (iatom < natoms); iatom++)
227         {
228             for (j = 0; (j < DIM); j++)
229             {
230                 row = DIM*iatom+j;
231                 for (k = 0; k < sparse_hessian->ndata[row]; k++)
232                 {
233                     col      = sparse_hessian->data[row][k].col;
234                     katom    = col/3;
235                     mass_fac = gmx_invsqrt(top->atoms.atom[iatom].m*top->atoms.atom[katom].m);
236                     sparse_hessian->data[row][k].value *= mass_fac;
237                 }
238             }
239         }
240     }
241     fprintf(stderr, "\nDiagonalizing to find eigenvectors 1 through %d...\n", neig);
242     fflush(stderr);
243
244     sparse_eigensolver(sparse_hessian, neig, eigenvalues, eigenvectors, 10000000);
245
246     /* Scale output eigenvectors */
247     if (bM && eigenvectors != NULL)
248     {
249         for (i = 0; i < neig; i++)
250         {
251             for (j = 0; j < natoms; j++)
252             {
253                 mass_fac = gmx_invsqrt(top->atoms.atom[j].m);
254                 for (k = 0; (k < DIM); k++)
255                 {
256                     eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
257                 }
258             }
259         }
260     }
261 }
262
263
264
265 int gmx_nmeig(int argc, char *argv[])
266 {
267     const char            *desc[] = {
268         "[THISMODULE] calculates the eigenvectors/values of a (Hessian) matrix,",
269         "which can be calculated with [gmx-mdrun].",
270         "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
271         "The structure is written first with t=0. The eigenvectors",
272         "are written as frames with the eigenvector number as timestamp.",
273         "The eigenvectors can be analyzed with [gmx-anaeig].",
274         "An ensemble of structures can be generated from the eigenvectors with",
275         "[gmx-nmens]. When mass weighting is used, the generated eigenvectors",
276         "will be scaled back to plain Cartesian coordinates before generating the",
277         "output. In this case, they will no longer be exactly orthogonal in the",
278         "standard Cartesian norm, but in the mass-weighted norm they would be.[PAR]",
279         "This program can be optionally used to compute quantum corrections to heat capacity",
280         "and enthalpy by providing an extra file argument [TT]-qcorr[tt]. See the GROMACS",
281         "manual, Chapter 1, for details. The result includes subtracting a harmonic",
282         "degree of freedom at the given temperature.",
283         "The total correction is printed on the terminal screen.",
284         "The recommended way of getting the corrections out is:[PAR]",
285         "[TT]gmx nmeig -s topol.tpr -f nm.mtx -first 7 -last 10000 -T 300 -qc [-constr][tt][PAR]",
286         "The [TT]-constr[tt] option should be used when bond constraints were used during the",
287         "simulation [BB]for all the covalent bonds[bb]. If this is not the case, ",
288         "you need to analyze the [TT]quant_corr.xvg[tt] file yourself.[PAR]",
289         "To make things more flexible, the program can also take virtual sites into account",
290         "when computing quantum corrections. When selecting [TT]-constr[tt] and",
291         "[TT]-qc[tt], the [TT]-begin[tt] and [TT]-end[tt] options will be set automatically as well.",
292         "Again, if you think you know it better, please check the [TT]eigenfreq.xvg[tt]",
293         "output."
294     };
295
296     static gmx_bool        bM    = TRUE, bCons = FALSE;
297     static int             begin = 1, end = 50, maxspec = 4000;
298     static real            T     = 298.15, width = 1;
299     t_pargs                pa[]  =
300     {
301         { "-m",  FALSE, etBOOL, {&bM},
302           "Divide elements of Hessian by product of sqrt(mass) of involved "
303           "atoms prior to diagonalization. This should be used for 'Normal Modes' "
304           "analysis" },
305         { "-first", FALSE, etINT, {&begin},
306           "First eigenvector to write away" },
307         { "-last",  FALSE, etINT, {&end},
308           "Last eigenvector to write away" },
309         { "-maxspec", FALSE, etINT, {&maxspec},
310           "Highest frequency (1/cm) to consider in the spectrum" },
311         { "-T",     FALSE, etREAL, {&T},
312           "Temperature for computing quantum heat capacity and enthalpy when using normal mode calculations to correct classical simulations" },
313         { "-constr", FALSE, etBOOL, {&bCons},
314           "If constraints were used in the simulation but not in the normal mode analysis (this is the recommended way of doing it) you will need to set this for computing the quantum corrections." },
315         { "-width",  FALSE, etREAL, {&width},
316           "Width (sigma) of the gaussian peaks (1/cm) when generating a spectrum" }
317     };
318     FILE                  *out, *qc, *spec;
319     int                    status, trjout;
320     t_topology             top;
321     gmx_mtop_t             mtop;
322     int                    ePBC;
323     rvec                  *top_x;
324     matrix                 box;
325     real                  *eigenvalues;
326     real                  *eigenvectors;
327     real                   rdum, mass_fac, qcvtot, qutot, qcv, qu;
328     int                    natoms, ndim, nrow, ncol, count, nharm, nvsite;
329     char                  *grpname;
330     int                    i, j, k, l, d, gnx;
331     gmx_bool               bSuck;
332     atom_id               *index;
333     t_tpxheader            tpx;
334     int                    version, generation;
335     real                   value, omega, nu;
336     real                   factor_gmx_to_omega2;
337     real                   factor_omega_to_wavenumber;
338     real                  *spectrum = NULL;
339     real                   wfac;
340     output_env_t           oenv;
341     const char            *qcleg[] = {
342         "Heat Capacity cV (J/mol K)",
343         "Enthalpy H (kJ/mol)"
344     };
345     real *                 full_hessian   = NULL;
346     gmx_sparsematrix_t *   sparse_hessian = NULL;
347
348     t_filenm               fnm[] = {
349         { efMTX, "-f", "hessian",    ffREAD  },
350         { efTPX, NULL, NULL,         ffREAD  },
351         { efXVG, "-of", "eigenfreq", ffWRITE },
352         { efXVG, "-ol", "eigenval",  ffWRITE },
353         { efXVG, "-os", "spectrum",  ffOPTWR },
354         { efXVG, "-qc", "quant_corr",  ffOPTWR },
355         { efTRN, "-v", "eigenvec",  ffWRITE }
356     };
357 #define NFILE asize(fnm)
358
359     if (!parse_common_args(&argc, argv, PCA_BE_NICE,
360                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
361     {
362         return 0;
363     }
364
365     /* Read tpr file for volume and number of harmonic terms */
366     read_tpxheader(ftp2fn(efTPX, NFILE, fnm), &tpx, TRUE, &version, &generation);
367     snew(top_x, tpx.natoms);
368
369     read_tpx(ftp2fn(efTPX, NFILE, fnm), NULL, box, &natoms,
370              top_x, NULL, NULL, &mtop);
371     if (bCons)
372     {
373         nharm = get_nharm(&mtop, &nvsite);
374     }
375     else
376     {
377         nharm  = 0;
378         nvsite = 0;
379     }
380     top = gmx_mtop_t_to_t_topology(&mtop);
381
382     bM   = TRUE;
383     ndim = DIM*natoms;
384
385     if (opt2bSet("-qc", NFILE, fnm))
386     {
387         begin = 7+DIM*nvsite;
388         end   = DIM*natoms;
389     }
390     if (begin < 1)
391     {
392         begin = 1;
393     }
394     if (end > ndim)
395     {
396         end = ndim;
397     }
398     printf("Using begin = %d and end = %d\n", begin, end);
399
400     /*open Hessian matrix */
401     gmx_mtxio_read(ftp2fn(efMTX, NFILE, fnm), &nrow, &ncol, &full_hessian, &sparse_hessian);
402
403     /* Memory for eigenvalues and eigenvectors (begin..end) */
404     snew(eigenvalues, nrow);
405     snew(eigenvectors, nrow*(end-begin+1));
406
407     /* If the Hessian is in sparse format we can calculate max (ndim-1) eigenvectors,
408      * and they must start at the first one. If this is not valid we convert to full matrix
409      * storage, but warn the user that we might run out of memory...
410      */
411     if ((sparse_hessian != NULL) && (begin != 1 || end == ndim))
412     {
413         if (begin != 1)
414         {
415             fprintf(stderr, "Cannot use sparse Hessian with first eigenvector != 1.\n");
416         }
417         else if (end == ndim)
418         {
419             fprintf(stderr, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
420         }
421
422         fprintf(stderr, "Will try to allocate memory and convert to full matrix representation...\n");
423
424         snew(full_hessian, nrow*ncol);
425         for (i = 0; i < nrow*ncol; i++)
426         {
427             full_hessian[i] = 0;
428         }
429
430         for (i = 0; i < sparse_hessian->nrow; i++)
431         {
432             for (j = 0; j < sparse_hessian->ndata[i]; j++)
433             {
434                 k                      = sparse_hessian->data[i][j].col;
435                 value                  = sparse_hessian->data[i][j].value;
436                 full_hessian[i*ndim+k] = value;
437                 full_hessian[k*ndim+i] = value;
438             }
439         }
440         gmx_sparsematrix_destroy(sparse_hessian);
441         sparse_hessian = NULL;
442         fprintf(stderr, "Converted sparse to full matrix storage.\n");
443     }
444
445     if (full_hessian != NULL)
446     {
447         /* Using full matrix storage */
448         nma_full_hessian(full_hessian, nrow, bM, &top, begin, end,
449                          eigenvalues, eigenvectors);
450     }
451     else
452     {
453         /* Sparse memory storage, allocate memory for eigenvectors */
454         snew(eigenvectors, ncol*end);
455         nma_sparse_hessian(sparse_hessian, bM, &top, end, eigenvalues, eigenvectors);
456     }
457
458     /* check the output, first 6 eigenvalues should be reasonably small */
459     bSuck = FALSE;
460     for (i = begin-1; (i < 6); i++)
461     {
462         if (fabs(eigenvalues[i]) > 1.0e-3)
463         {
464             bSuck = TRUE;
465         }
466     }
467     if (bSuck)
468     {
469         fprintf(stderr, "\nOne of the lowest 6 eigenvalues has a non-zero value.\n");
470         fprintf(stderr, "This could mean that the reference structure was not\n");
471         fprintf(stderr, "properly energy minimized.\n");
472     }
473
474     /* now write the output */
475     fprintf (stderr, "Writing eigenvalues...\n");
476     out = xvgropen(opt2fn("-ol", NFILE, fnm),
477                    "Eigenvalues", "Eigenvalue index", "Eigenvalue [Gromacs units]",
478                    oenv);
479     if (output_env_get_print_xvgr_codes(oenv))
480     {
481         if (bM)
482         {
483             fprintf(out, "@ subtitle \"mass weighted\"\n");
484         }
485         else
486         {
487             fprintf(out, "@ subtitle \"not mass weighted\"\n");
488         }
489     }
490
491     for (i = 0; i <= (end-begin); i++)
492     {
493         fprintf (out, "%6d %15g\n", begin+i, eigenvalues[i]);
494     }
495     gmx_ffclose(out);
496
497
498     if (opt2bSet("-qc", NFILE, fnm))
499     {
500         qc = xvgropen(opt2fn("-qc", NFILE, fnm), "Quantum Corrections", "Eigenvector index", "", oenv);
501         xvgr_legend(qc, asize(qcleg), qcleg, oenv);
502         qcvtot = qutot = 0;
503     }
504     else
505     {
506         qc = NULL;
507     }
508     printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
509
510     out = xvgropen(opt2fn("-of", NFILE, fnm),
511                    "Eigenfrequencies", "Eigenvector index", "Wavenumber [cm\\S-1\\N]",
512                    oenv);
513     if (output_env_get_print_xvgr_codes(oenv))
514     {
515         if (bM)
516         {
517             fprintf(out, "@ subtitle \"mass weighted\"\n");
518         }
519         else
520         {
521             fprintf(out, "@ subtitle \"not mass weighted\"\n");
522         }
523     }
524     /* Spectrum ? */
525     spec = NULL;
526     if (opt2bSet("-os", NFILE, fnm) && (maxspec > 0))
527     {
528         snew(spectrum, maxspec);
529         spec = xvgropen(opt2fn("-os", NFILE, fnm),
530                         "Vibrational spectrum based on harmonic approximation",
531                         "\\f{12}w\\f{4} (cm\\S-1\\N)",
532                         "Intensity [Gromacs units]",
533                         oenv);
534         for (i = 0; (i < maxspec); i++)
535         {
536             spectrum[i] = 0;
537         }
538     }
539
540     /* Gromacs units are kJ/(mol*nm*nm*amu),
541      * where amu is the atomic mass unit.
542      *
543      * For the eigenfrequencies we want to convert this to spectroscopic absorption
544      * wavenumbers given in cm^(-1), which is the frequency divided by the speed of
545      * light. Do this by first converting to omega^2 (units 1/s), take the square
546      * root, and finally divide by the speed of light (nm/ps in gromacs).
547      */
548     factor_gmx_to_omega2       = 1.0E21/(AVOGADRO*AMU);
549     factor_omega_to_wavenumber = 1.0E-5/(2.0*M_PI*SPEED_OF_LIGHT);
550
551     for (i = begin; (i <= end); i++)
552     {
553         value = eigenvalues[i-begin];
554         if (value < 0)
555         {
556             value = 0;
557         }
558         omega = sqrt(value*factor_gmx_to_omega2);
559         nu    = 1e-12*omega/(2*M_PI);
560         value = omega*factor_omega_to_wavenumber;
561         fprintf (out, "%6d %15g\n", i, value);
562         if (NULL != spec)
563         {
564             wfac = eigenvalues[i-begin]/(width*sqrt(2*M_PI));
565             for (j = 0; (j < maxspec); j++)
566             {
567                 spectrum[j] += wfac*exp(-sqr(j-value)/(2*sqr(width)));
568             }
569         }
570         if (NULL != qc)
571         {
572             qcv = cv_corr(nu, T);
573             qu  = u_corr(nu, T);
574             if (i > end-nharm)
575             {
576                 qcv += BOLTZ*KILO;
577                 qu  += BOLTZ*T;
578             }
579             fprintf (qc, "%6d %15g %15g\n", i, qcv, qu);
580             qcvtot += qcv;
581             qutot  += qu;
582         }
583     }
584     gmx_ffclose(out);
585     if (NULL != spec)
586     {
587         for (j = 0; (j < maxspec); j++)
588         {
589             fprintf(spec, "%10g  %10g\n", 1.0*j, spectrum[j]);
590         }
591         gmx_ffclose(spec);
592     }
593     if (NULL != qc)
594     {
595         printf("Quantum corrections for harmonic degrees of freedom\n");
596         printf("Use appropriate -first and -last options to get reliable results.\n");
597         printf("There were %d constraints and %d vsites in the simulation\n",
598                nharm, nvsite);
599         printf("Total correction to cV = %g J/mol K\n", qcvtot);
600         printf("Total correction to  H = %g kJ/mol\n", qutot);
601         gmx_ffclose(qc);
602         please_cite(stdout, "Caleman2011b");
603     }
604     /* Writing eigenvectors. Note that if mass scaling was used, the eigenvectors
605      * were scaled back from mass weighted cartesian to plain cartesian in the
606      * nma_full_hessian() or nma_sparse_hessian() routines. Mass scaled vectors
607      * will not be strictly orthogonal in plain cartesian scalar products.
608      */
609     write_eigenvectors(opt2fn("-v", NFILE, fnm), natoms, eigenvectors, FALSE, begin, end,
610                        eWXR_NO, NULL, FALSE, top_x, bM, eigenvalues);
611
612     return 0;
613 }