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