3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
58 #include "mtop_util.h"
63 #include "gromacs/linearalgebra/eigensolver.h"
64 #include "gromacs/linearalgebra/mtxio.h"
65 #include "gromacs/linearalgebra/sparsematrix.h"
67 static double cv_corr(double nu, double T)
69 double x = PLANCK*nu/(BOLTZ*T);
78 return BOLTZ*KILO*(ex*sqr(x)/sqr(ex-1) - 1);
82 static double u_corr(double nu, double T)
84 double x = PLANCK*nu/(BOLTZ*T);
93 return BOLTZ*T*(0.5*x - 1 + x/(ex-1));
97 static int get_nharm_mt(gmx_moltype_t *mt)
99 static int harm_func[] = { F_BONDS };
103 for (i = 0; (i < asize(harm_func)); i++)
106 nh += mt->ilist[ft].nr/(interaction_function[ft].nratoms+1);
111 static int get_nvsite_mt(gmx_moltype_t *mt)
113 static int vs_func[] = {
114 F_VSITE2, F_VSITE3, F_VSITE3FD, F_VSITE3FAD,
115 F_VSITE3OUT, F_VSITE4FD, F_VSITE4FDN, F_VSITEN
120 for (i = 0; (i < asize(vs_func)); i++)
123 nh += mt->ilist[ft].nr/(interaction_function[ft].nratoms+1);
128 static int get_nharm(gmx_mtop_t *mtop, int *nvsites)
134 for (j = 0; (j < mtop->nmolblock); j++)
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]));
145 nma_full_hessian(real * hess,
158 natoms = top->atoms.nr;
160 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
164 for (i = 0; (i < natoms); i++)
166 for (j = 0; (j < DIM); j++)
168 for (k = 0; (k < natoms); k++)
170 mass_fac = gmx_invsqrt(top->atoms.atom[i].m*top->atoms.atom[k].m);
171 for (l = 0; (l < DIM); l++)
173 hess[(i*DIM+j)*ndim+k*DIM+l] *= mass_fac;
180 /* call diagonalization routine. */
182 fprintf(stderr, "\nDiagonalizing to find vectors %d through %d...\n", begin, end);
185 eigensolver(hess, ndim, begin-1, end-1, eigenvalues, eigenvectors);
187 /* And scale the output eigenvectors */
188 if (bM && eigenvectors != NULL)
190 for (i = 0; i < (end-begin+1); i++)
192 for (j = 0; j < natoms; j++)
194 mass_fac = gmx_invsqrt(top->atoms.atom[j].m);
195 for (k = 0; (k < DIM); k++)
197 eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
207 nma_sparse_hessian(gmx_sparsematrix_t * sparse_hessian,
221 natoms = top->atoms.nr;
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 */
229 for (iatom = 0; (iatom < natoms); iatom++)
231 for (j = 0; (j < DIM); j++)
234 for (k = 0; k < sparse_hessian->ndata[row]; k++)
236 col = sparse_hessian->data[row][k].col;
238 mass_fac = gmx_invsqrt(top->atoms.atom[iatom].m*top->atoms.atom[katom].m);
239 sparse_hessian->data[row][k].value *= mass_fac;
244 fprintf(stderr, "\nDiagonalizing to find eigenvectors 1 through %d...\n", neig);
247 sparse_eigensolver(sparse_hessian, neig, eigenvalues, eigenvectors, 10000000);
249 /* Scale output eigenvectors */
250 if (bM && eigenvectors != NULL)
252 for (i = 0; i < neig; i++)
254 for (j = 0; j < natoms; j++)
256 mass_fac = gmx_invsqrt(top->atoms.atom[j].m);
257 for (k = 0; (k < DIM); k++)
259 eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
268 int gmx_nmeig(int argc, char *argv[])
270 const char *desc[] = {
271 "[TT]g_nmeig[tt] calculates the eigenvectors/values of a (Hessian) matrix,",
272 "which can be calculated with [TT]mdrun[tt].",
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 [TT]g_anaeig[tt].",
277 "An ensemble of structures can be generated from the eigenvectors with",
278 "[TT]g_nmens[tt]. 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]g_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]",
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;
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' "
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" }
321 FILE *out, *qc, *spec;
330 real rdum, mass_fac, qcvtot, qutot, qcv, qu;
331 int natoms, ndim, nrow, ncol, count, nharm, nvsite;
333 int i, j, k, l, d, gnx;
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;
344 const char *qcleg[] = {
345 "Heat Capacity cV (J/mol K)",
346 "Enthalpy H (kJ/mol)"
348 real * full_hessian = NULL;
349 gmx_sparsematrix_t * sparse_hessian = NULL;
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 }
360 #define NFILE asize(fnm)
362 parse_common_args(&argc, argv, PCA_BE_NICE,
363 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
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);
369 read_tpx(ftp2fn(efTPX, NFILE, fnm), NULL, box, &natoms,
370 top_x, NULL, NULL, &mtop);
373 nharm = get_nharm(&mtop, &nvsite);
380 top = gmx_mtop_t_to_t_topology(&mtop);
385 if (opt2bSet("-qc", NFILE, fnm))
387 begin = 7+DIM*nvsite;
398 printf("Using begin = %d and end = %d\n", begin, end);
400 /*open Hessian matrix */
401 gmx_mtxio_read(ftp2fn(efMTX, NFILE, fnm), &nrow, &ncol, &full_hessian, &sparse_hessian);
403 /* Memory for eigenvalues and eigenvectors (begin..end) */
404 snew(eigenvalues, nrow);
405 snew(eigenvectors, nrow*(end-begin+1));
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...
411 if ((sparse_hessian != NULL) && (begin != 1 || end == ndim))
415 fprintf(stderr, "Cannot use sparse Hessian with first eigenvector != 1.\n");
417 else if (end == ndim)
419 fprintf(stderr, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
422 fprintf(stderr, "Will try to allocate memory and convert to full matrix representation...\n");
424 snew(full_hessian, nrow*ncol);
425 for (i = 0; i < nrow*ncol; i++)
430 for (i = 0; i < sparse_hessian->nrow; i++)
432 for (j = 0; j < sparse_hessian->ndata[i]; j++)
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;
440 gmx_sparsematrix_destroy(sparse_hessian);
441 sparse_hessian = NULL;
442 fprintf(stderr, "Converted sparse to full matrix storage.\n");
445 if (full_hessian != NULL)
447 /* Using full matrix storage */
448 nma_full_hessian(full_hessian, nrow, bM, &top, begin, end,
449 eigenvalues, eigenvectors);
453 /* Sparse memory storage, allocate memory for eigenvectors */
454 snew(eigenvectors, ncol*end);
455 nma_sparse_hessian(sparse_hessian, bM, &top, end, eigenvalues, eigenvectors);
458 /* check the output, first 6 eigenvalues should be reasonably small */
460 for (i = begin-1; (i < 6); i++)
462 if (fabs(eigenvalues[i]) > 1.0e-3)
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");
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]",
479 if (output_env_get_print_xvgr_codes(oenv))
483 fprintf(out, "@ subtitle \"mass weighted\"\n");
487 fprintf(out, "@ subtitle \"not mass weighted\"\n");
491 for (i = 0; i <= (end-begin); i++)
493 fprintf (out, "%6d %15g\n", begin+i, eigenvalues[i]);
498 if (opt2bSet("-qc", NFILE, fnm))
500 qc = xvgropen(opt2fn("-qc", NFILE, fnm), "Quantum Corrections", "Eigenvector index", "", oenv);
501 xvgr_legend(qc, asize(qcleg), qcleg, oenv);
508 printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
510 out = xvgropen(opt2fn("-of", NFILE, fnm),
511 "Eigenfrequencies", "Eigenvector index", "Wavenumber [cm\\S-1\\N]",
513 if (output_env_get_print_xvgr_codes(oenv))
517 fprintf(out, "@ subtitle \"mass weighted\"\n");
521 fprintf(out, "@ subtitle \"not mass weighted\"\n");
526 if (opt2bSet("-os", NFILE, fnm) && (maxspec > 0))
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]",
534 for (i = 0; (i < maxspec); i++)
540 /* Gromacs units are kJ/(mol*nm*nm*amu),
541 * where amu is the atomic mass unit.
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).
548 factor_gmx_to_omega2 = 1.0E21/(AVOGADRO*AMU);
549 factor_omega_to_wavenumber = 1.0E-5/(2.0*M_PI*SPEED_OF_LIGHT);
551 for (i = begin; (i <= end); i++)
553 value = eigenvalues[i-begin];
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);
564 wfac = eigenvalues[i-begin]/(width*sqrt(2*M_PI));
565 for (j = 0; (j < maxspec); j++)
567 spectrum[j] += wfac*exp(-sqr(j-value)/(2*sqr(width)));
572 qcv = cv_corr(nu, T);
579 fprintf (qc, "%6d %15g %15g\n", i, qcv, qu);
587 for (j = 0; (j < maxspec); j++)
589 fprintf(spec, "%10g %10g\n", 1.0*j, spectrum[j]);
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",
599 printf("Total correction to cV = %g J/mol K\n", qcvtot);
600 printf("Total correction to H = %g kJ/mol\n", qutot);
602 please_cite(stdout, "Caleman2011b");
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.
609 write_eigenvectors(opt2fn("-v", NFILE, fnm), natoms, eigenvectors, FALSE, begin, end,
610 eWXR_NO, NULL, FALSE, top_x, bM, eigenvalues);