2 * This file is part of the GROMACS molecular simulation package.
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,2016,2017,2018, 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.
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.
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.
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.
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.
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.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/fileio/mtxio.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/eigio.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/gstat.h"
52 #include "gromacs/linearalgebra/eigensolver.h"
53 #include "gromacs/linearalgebra/sparsematrix.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/mtop_util.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/pleasecite.h"
65 #include "gromacs/utility/smalloc.h"
67 #include "thermochemistry.h"
69 static double cv_corr(double nu, double T)
71 double x = PLANCK*nu/(BOLTZ*T);
72 double ex = std::exp(x);
80 return BOLTZ*KILO*(ex*gmx::square(x)/gmx::square(ex-1) - 1);
84 static double u_corr(double nu, double T)
86 double x = PLANCK*nu/(BOLTZ*T);
87 double ex = std::exp(x);
95 return BOLTZ*T*(0.5*x - 1 + x/(ex-1));
99 static size_t get_nharm_mt(const gmx_moltype_t *mt)
101 static int harm_func[] = { F_BONDS };
105 for (i = 0; (i < asize(harm_func)); i++)
108 nh += mt->ilist[ft].nr/(interaction_function[ft].nratoms+1);
113 static int get_nharm(const gmx_mtop_t *mtop)
117 for (const gmx_molblock_t &molb : mtop->molblock)
119 nh += molb.nmol * get_nharm_mt(&(mtop->moltype[molb.type]));
125 nma_full_hessian(real *hess,
128 const t_topology *top,
129 const std::vector<size_t> &atom_index,
137 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
141 for (size_t i = 0; (i < atom_index.size()); i++)
143 size_t ai = atom_index[i];
144 for (size_t j = 0; (j < DIM); j++)
146 for (size_t k = 0; (k < atom_index.size()); k++)
148 size_t ak = atom_index[k];
149 mass_fac = gmx::invsqrt(top->atoms.atom[ai].m*top->atoms.atom[ak].m);
150 for (size_t l = 0; (l < DIM); l++)
152 hess[(i*DIM+j)*ndim+k*DIM+l] *= mass_fac;
159 /* call diagonalization routine. */
161 fprintf(stderr, "\nDiagonalizing to find vectors %d through %d...\n", begin, end);
164 eigensolver(hess, ndim, begin-1, end-1, eigenvalues, eigenvectors);
166 /* And scale the output eigenvectors */
167 if (bM && eigenvectors != nullptr)
169 for (int i = 0; i < (end-begin+1); i++)
171 for (size_t j = 0; j < atom_index.size(); j++)
173 size_t aj = atom_index[j];
174 mass_fac = gmx::invsqrt(top->atoms.atom[aj].m);
175 for (size_t k = 0; (k < DIM); k++)
177 eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
187 nma_sparse_hessian(gmx_sparsematrix_t *sparse_hessian,
189 const t_topology *top,
190 const std::vector<size_t> &atom_index,
201 ndim = DIM*atom_index.size();
203 /* Cannot check symmetry since we only store half matrix */
204 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
206 GMX_RELEASE_ASSERT(sparse_hessian != nullptr, "NULL matrix pointer provided to nma_sparse_hessian");
210 for (size_t iatom = 0; (iatom < atom_index.size()); iatom++)
212 size_t ai = atom_index[iatom];
213 for (size_t j = 0; (j < DIM); j++)
216 for (k = 0; k < sparse_hessian->ndata[row]; k++)
218 col = sparse_hessian->data[row][k].col;
220 size_t ak = atom_index[katom];
221 mass_fac = gmx::invsqrt(top->atoms.atom[ai].m*top->atoms.atom[ak].m);
222 sparse_hessian->data[row][k].value *= mass_fac;
227 fprintf(stderr, "\nDiagonalizing to find eigenvectors 1 through %d...\n", neig);
230 sparse_eigensolver(sparse_hessian, neig, eigenvalues, eigenvectors, 10000000);
232 /* Scale output eigenvectors */
233 if (bM && eigenvectors != nullptr)
235 for (i = 0; i < neig; i++)
237 for (size_t j = 0; j < atom_index.size(); j++)
239 size_t aj = atom_index[j];
240 mass_fac = gmx::invsqrt(top->atoms.atom[aj].m);
241 for (k = 0; (k < DIM); k++)
243 eigenvectors[i*ndim+j*DIM+k] *= mass_fac;
251 /* Returns a pointer for eigenvector storage */
252 static real *allocateEigenvectors(int nrow, int first, int last,
262 numVector = last - first + 1;
264 size_t vectorsSize = static_cast<size_t>(nrow)*static_cast<size_t>(numVector);
265 /* We can't have more than INT_MAX elements.
266 * Relaxing this restriction probably requires changing lots of loop
267 * variable types in the linear algebra code.
269 if (vectorsSize > INT_MAX)
271 gmx_fatal(FARGS, "You asked to store %d eigenvectors of size %d, which requires more than the supported %d elements; %sdecrease -last",
272 numVector, nrow, INT_MAX,
273 ignoreBegin ? "" : "increase -first and/or ");
277 snew(eigenvectors, vectorsSize);
283 int gmx_nmeig(int argc, char *argv[])
285 const char *desc[] = {
286 "[THISMODULE] calculates the eigenvectors/values of a (Hessian) matrix,",
287 "which can be calculated with [gmx-mdrun].",
288 "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
289 "The structure is written first with t=0. The eigenvectors",
290 "are written as frames with the eigenvector number and eigenvalue",
291 "written as step number and timestamp, respectively.",
292 "The eigenvectors can be analyzed with [gmx-anaeig].",
293 "An ensemble of structures can be generated from the eigenvectors with",
294 "[gmx-nmens]. When mass weighting is used, the generated eigenvectors",
295 "will be scaled back to plain Cartesian coordinates before generating the",
296 "output. In this case, they will no longer be exactly orthogonal in the",
297 "standard Cartesian norm, but in the mass-weighted norm they would be.[PAR]",
298 "This program can be optionally used to compute quantum corrections to heat capacity",
299 "and enthalpy by providing an extra file argument [TT]-qcorr[tt]. See the GROMACS",
300 "manual, Chapter 1, for details. The result includes subtracting a harmonic",
301 "degree of freedom at the given temperature.",
302 "The total correction is printed on the terminal screen.",
303 "The recommended way of getting the corrections out is:[PAR]",
304 "[TT]gmx nmeig -s topol.tpr -f nm.mtx -first 7 -last 10000 -T 300 -qc [-constr][tt][PAR]",
305 "The [TT]-constr[tt] option should be used when bond constraints were used during the",
306 "simulation [BB]for all the covalent bonds[bb]. If this is not the case, ",
307 "you need to analyze the [TT]quant_corr.xvg[tt] file yourself.[PAR]",
308 "To make things more flexible, the program can also take virtual sites into account",
309 "when computing quantum corrections. When selecting [TT]-constr[tt] and",
310 "[TT]-qc[tt], the [TT]-begin[tt] and [TT]-end[tt] options will be set automatically as well.",
311 "Again, if you think you know it better, please check the [TT]eigenfreq.xvg[tt]",
315 static gmx_bool bM = TRUE, bCons = FALSE, bLinear = FALSE;
316 static int begin = 1, end = 50, maxspec = 4000;
317 static real T = 298.15, width = 1;
320 { "-m", FALSE, etBOOL, {&bM},
321 "Divide elements of Hessian by product of sqrt(mass) of involved "
322 "atoms prior to diagonalization. This should be used for 'Normal Modes' "
324 { "-linear", FALSE, etBOOL, {&bLinear},
325 "This should be set in order to get correct entropies for linear molecules" },
326 { "-first", FALSE, etINT, {&begin},
327 "First eigenvector to write away" },
328 { "-last", FALSE, etINT, {&end},
329 "Last eigenvector to write away. -1 (default) is use all dimensions." },
330 { "-maxspec", FALSE, etINT, {&maxspec},
331 "Highest frequency (1/cm) to consider in the spectrum" },
332 { "-T", FALSE, etREAL, {&T},
333 "Temperature for computing quantum heat capacity and enthalpy when using normal mode calculations to correct classical simulations" },
334 { "-constr", FALSE, etBOOL, {&bCons},
335 "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." },
336 { "-width", FALSE, etREAL, {&width},
337 "Width (sigma) of the gaussian peaks (1/cm) when generating a spectrum" }
339 FILE *out, *qc, *spec;
346 real qcvtot, qutot, qcv, qu;
349 real value, omega, nu;
350 real factor_gmx_to_omega2;
351 real factor_omega_to_wavenumber;
352 real *spectrum = nullptr;
354 gmx_output_env_t *oenv;
355 const char *qcleg[] = {
356 "Heat Capacity cV (J/mol K)",
357 "Enthalpy H (kJ/mol)"
359 real * full_hessian = nullptr;
360 gmx_sparsematrix_t * sparse_hessian = nullptr;
363 { efMTX, "-f", "hessian", ffREAD },
364 { efTPR, nullptr, nullptr, ffREAD },
365 { efXVG, "-of", "eigenfreq", ffWRITE },
366 { efXVG, "-ol", "eigenval", ffWRITE },
367 { efXVG, "-os", "spectrum", ffOPTWR },
368 { efXVG, "-qc", "quant_corr", ffOPTWR },
369 { efTRN, "-v", "eigenvec", ffWRITE }
371 #define NFILE asize(fnm)
373 if (!parse_common_args(&argc, argv, 0,
374 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
379 /* Read tpr file for volume and number of harmonic terms */
380 read_tpxheader(ftp2fn(efTPR, NFILE, fnm), &tpx, TRUE);
381 snew(top_x, tpx.natoms);
384 read_tpx(ftp2fn(efTPR, NFILE, fnm), nullptr, box, &natoms_tpx,
385 top_x, nullptr, &mtop);
389 nharm = get_nharm(&mtop);
391 std::vector<size_t> atom_index = get_atom_index(&mtop);
393 top = gmx_mtop_t_to_t_topology(&mtop, true);
396 int ndim = DIM*atom_index.size();
398 if (opt2bSet("-qc", NFILE, fnm))
407 if (end == -1 || end > ndim)
411 printf("Using begin = %d and end = %d\n", begin, end);
413 /*open Hessian matrix */
415 gmx_mtxio_read(ftp2fn(efMTX, NFILE, fnm), &nrow, &ncol, &full_hessian, &sparse_hessian);
417 /* If the Hessian is in sparse format we can calculate max (ndim-1) eigenvectors,
418 * If this is not valid we convert to full matrix storage,
419 * but warn the user that we might run out of memory...
421 if ((sparse_hessian != nullptr) && (end == ndim))
423 fprintf(stderr, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
425 fprintf(stderr, "Will try to allocate memory and convert to full matrix representation...\n");
427 size_t hessianSize = static_cast<size_t>(nrow)*static_cast<size_t>(ncol);
428 /* Allowing Hessians larger than INT_MAX probably only makes sense
429 * with (OpenMP) parallel diagonalization routines, since with a single
430 * thread it will takes months.
432 if (hessianSize > INT_MAX)
434 gmx_fatal(FARGS, "Hessian size is %d x %d, which is larger than the maximum allowed %d elements.",
435 nrow, ncol, INT_MAX);
437 snew(full_hessian, hessianSize);
438 for (i = 0; i < nrow*ncol; i++)
443 for (i = 0; i < sparse_hessian->nrow; i++)
445 for (j = 0; j < sparse_hessian->ndata[i]; j++)
447 k = sparse_hessian->data[i][j].col;
448 value = sparse_hessian->data[i][j].value;
449 full_hessian[i*ndim+k] = value;
450 full_hessian[k*ndim+i] = value;
453 gmx_sparsematrix_destroy(sparse_hessian);
454 sparse_hessian = nullptr;
455 fprintf(stderr, "Converted sparse to full matrix storage.\n");
458 snew(eigenvalues, nrow);
460 if (full_hessian != nullptr)
462 /* Using full matrix storage */
463 eigenvectors = allocateEigenvectors(nrow, begin, end, false);
465 nma_full_hessian(full_hessian, nrow, bM, &top, atom_index, begin, end,
466 eigenvalues, eigenvectors);
470 assert(sparse_hessian);
471 /* Sparse memory storage, allocate memory for eigenvectors */
472 eigenvectors = allocateEigenvectors(nrow, begin, end, true);
474 nma_sparse_hessian(sparse_hessian, bM, &top, atom_index, end, eigenvalues, eigenvectors);
477 /* check the output, first 6 eigenvalues should be reasonably small */
478 gmx_bool bSuck = FALSE;
479 for (i = begin-1; (i < 6); i++)
481 if (std::abs(eigenvalues[i]) > 1.0e-3)
488 fprintf(stderr, "\nOne of the lowest 6 eigenvalues has a non-zero value.\n");
489 fprintf(stderr, "This could mean that the reference structure was not\n");
490 fprintf(stderr, "properly energy minimized.\n");
493 /* now write the output */
494 fprintf (stderr, "Writing eigenvalues...\n");
495 out = xvgropen(opt2fn("-ol", NFILE, fnm),
496 "Eigenvalues", "Eigenvalue index", "Eigenvalue [Gromacs units]",
498 if (output_env_get_print_xvgr_codes(oenv))
502 fprintf(out, "@ subtitle \"mass weighted\"\n");
506 fprintf(out, "@ subtitle \"not mass weighted\"\n");
510 for (i = 0; i <= (end-begin); i++)
512 fprintf (out, "%6d %15g\n", begin+i, eigenvalues[i]);
517 if (opt2bSet("-qc", NFILE, fnm))
519 qc = xvgropen(opt2fn("-qc", NFILE, fnm), "Quantum Corrections", "Eigenvector index", "", oenv);
520 xvgr_legend(qc, asize(qcleg), qcleg, oenv);
527 printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
529 out = xvgropen(opt2fn("-of", NFILE, fnm),
530 "Eigenfrequencies", "Eigenvector index", "Wavenumber [cm\\S-1\\N]",
532 if (output_env_get_print_xvgr_codes(oenv))
536 fprintf(out, "@ subtitle \"mass weighted\"\n");
540 fprintf(out, "@ subtitle \"not mass weighted\"\n");
545 if (opt2bSet("-os", NFILE, fnm) && (maxspec > 0))
547 snew(spectrum, maxspec);
548 spec = xvgropen(opt2fn("-os", NFILE, fnm),
549 "Vibrational spectrum based on harmonic approximation",
550 "\\f{12}w\\f{4} (cm\\S-1\\N)",
551 "Intensity [Gromacs units]",
553 for (i = 0; (i < maxspec); i++)
559 /* Gromacs units are kJ/(mol*nm*nm*amu),
560 * where amu is the atomic mass unit.
562 * For the eigenfrequencies we want to convert this to spectroscopic absorption
563 * wavenumbers given in cm^(-1), which is the frequency divided by the speed of
564 * light. Do this by first converting to omega^2 (units 1/s), take the square
565 * root, and finally divide by the speed of light (nm/ps in gromacs).
567 factor_gmx_to_omega2 = 1.0E21/(AVOGADRO*AMU);
568 factor_omega_to_wavenumber = 1.0E-5/(2.0*M_PI*SPEED_OF_LIGHT);
570 for (i = begin; (i <= end); i++)
572 value = eigenvalues[i-begin];
577 omega = std::sqrt(value*factor_gmx_to_omega2);
578 nu = 1e-12*omega/(2*M_PI);
579 value = omega*factor_omega_to_wavenumber;
580 fprintf (out, "%6d %15g\n", i, value);
583 wfac = eigenvalues[i-begin]/(width*std::sqrt(2*M_PI));
584 for (j = 0; (j < maxspec); j++)
586 spectrum[j] += wfac*std::exp(-gmx::square(j-value)/(2*gmx::square(width)));
591 qcv = cv_corr(nu, T);
598 fprintf (qc, "%6d %15g %15g\n", i, qcv, qu);
606 for (j = 0; (j < maxspec); j++)
608 fprintf(spec, "%10g %10g\n", 1.0*j, spectrum[j]);
614 printf("Quantum corrections for harmonic degrees of freedom\n");
615 printf("Use appropriate -first and -last options to get reliable results.\n");
616 printf("There were %d constraints in the simulation\n", nharm);
617 printf("Total correction to cV = %g J/mol K\n", qcvtot);
618 printf("Total correction to H = %g kJ/mol\n", qutot);
620 please_cite(stdout, "Caleman2011b");
622 /* Writing eigenvectors. Note that if mass scaling was used, the eigenvectors
623 * were scaled back from mass weighted cartesian to plain cartesian in the
624 * nma_full_hessian() or nma_sparse_hessian() routines. Mass scaled vectors
625 * will not be strictly orthogonal in plain cartesian scalar products.
627 const real *eigenvectorPtr;
628 if (full_hessian != nullptr)
630 eigenvectorPtr = eigenvectors;
634 /* The sparse matrix diagonalization store all eigenvectors up to end */
635 eigenvectorPtr = eigenvectors + (begin - 1)*atom_index.size();
637 write_eigenvectors(opt2fn("-v", NFILE, fnm), atom_index.size(), eigenvectorPtr, FALSE, begin, end,
638 eWXR_NO, nullptr, FALSE, top_x, bM, eigenvalues);
642 printf("The Entropy due to the Quasi Harmonic approximation is %g J/mol K\n",
643 calc_entropy_quasi_harmonic(DIM*atom_index.size(),
644 eigenvalues, T, bLinear));
648 printf("Cannot compute entropy when -first = %d\n", begin);