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,2019, 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/gmxana/princ.h"
53 #include "gromacs/linearalgebra/eigensolver.h"
54 #include "gromacs/linearalgebra/sparsematrix.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/topology/ifunc.h"
60 #include "gromacs/topology/mtop_util.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/arraysize.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/futil.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/pleasecite.h"
67 #include "gromacs/utility/smalloc.h"
69 #include "thermochemistry.h"
71 static double cv_corr(double nu, double T)
73 double x = PLANCK * nu / (BOLTZ * T);
74 double ex = std::exp(x);
82 return BOLTZ * KILO * (ex * gmx::square(x) / gmx::square(ex - 1) - 1);
86 static double u_corr(double nu, double T)
88 double x = PLANCK * nu / (BOLTZ * T);
89 double ex = std::exp(x);
97 return BOLTZ * T * (0.5 * x - 1 + x / (ex - 1));
101 static size_t get_nharm_mt(const gmx_moltype_t* mt)
103 static int harm_func[] = { F_BONDS };
107 for (i = 0; (i < asize(harm_func)); i++)
110 nh += mt->ilist[ft].size() / (interaction_function[ft].nratoms + 1);
115 static int get_nharm(const gmx_mtop_t* mtop)
119 for (const gmx_molblock_t& molb : mtop->molblock)
121 nh += molb.nmol * get_nharm_mt(&(mtop->moltype[molb.type]));
126 static void nma_full_hessian(real* hess,
129 const t_topology* top,
130 gmx::ArrayRef<const int> atom_index,
138 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
142 for (int i = 0; (i < atom_index.ssize()); i++)
144 size_t ai = atom_index[i];
145 for (size_t j = 0; (j < DIM); j++)
147 for (int k = 0; (k < atom_index.ssize()); k++)
149 size_t ak = atom_index[k];
150 mass_fac = gmx::invsqrt(top->atoms.atom[ai].m * top->atoms.atom[ak].m);
151 for (size_t l = 0; (l < DIM); l++)
153 hess[(i * DIM + j) * ndim + k * DIM + l] *= mass_fac;
160 /* call diagonalization routine. */
162 fprintf(stderr, "\nDiagonalizing to find vectors %d through %d...\n", begin, end);
165 eigensolver(hess, ndim, begin - 1, end - 1, eigenvalues, eigenvectors);
167 /* And scale the output eigenvectors */
168 if (bM && eigenvectors != nullptr)
170 for (int i = 0; i < (end - begin + 1); i++)
172 for (gmx::index j = 0; j < atom_index.ssize(); j++)
174 size_t aj = atom_index[j];
175 mass_fac = gmx::invsqrt(top->atoms.atom[aj].m);
176 for (size_t k = 0; (k < DIM); k++)
178 eigenvectors[i * ndim + j * DIM + k] *= mass_fac;
186 static void nma_sparse_hessian(gmx_sparsematrix_t* sparse_hessian,
188 const t_topology* top,
189 gmx::ArrayRef<const int> atom_index,
200 ndim = DIM * atom_index.size();
202 /* Cannot check symmetry since we only store half matrix */
203 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
205 GMX_RELEASE_ASSERT(sparse_hessian != nullptr,
206 "NULL matrix pointer provided to nma_sparse_hessian");
210 for (int iatom = 0; (iatom < atom_index.ssize()); iatom++)
212 size_t ai = atom_index[iatom];
213 for (size_t j = 0; (j < DIM); j++)
215 row = DIM * iatom + 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 (gmx::index j = 0; j < atom_index.ssize(); 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, bool ignoreBegin)
261 numVector = last - first + 1;
263 size_t vectorsSize = static_cast<size_t>(nrow) * static_cast<size_t>(numVector);
264 /* We can't have more than INT_MAX elements.
265 * Relaxing this restriction probably requires changing lots of loop
266 * variable types in the linear algebra code.
268 if (vectorsSize > INT_MAX)
271 "You asked to store %d eigenvectors of size %d, which requires more than the "
272 "supported %d elements; %sdecrease -last",
273 numVector, nrow, INT_MAX, ignoreBegin ? "" : "increase -first and/or ");
277 snew(eigenvectors, vectorsSize);
282 /*! \brief Compute heat capacity due to translational motion
284 static double calcTranslationalHeatCapacity()
289 /*! \brief Compute internal energy due to translational motion
291 static double calcTranslationalInternalEnergy(double T)
293 return BOLTZ * T * 1.5;
296 /*! \brief Compute heat capacity due to rotational motion
298 * \param[in] linear Should be TRUE if this is a linear molecule
299 * \param[in] T Temperature
300 * \return The heat capacity at constant volume
302 static double calcRotationalInternalEnergy(gmx_bool linear, double T)
310 return BOLTZ * T * 1.5;
314 /*! \brief Compute heat capacity due to rotational motion
316 * \param[in] linear Should be TRUE if this is a linear molecule
317 * \return The heat capacity at constant volume
319 static double calcRotationalHeatCapacity(gmx_bool linear)
331 static void analyzeThermochemistry(FILE* fp,
332 const t_topology& top,
334 gmx::ArrayRef<const int> atom_index,
342 std::vector<int> index(atom_index.begin(), atom_index.end());
345 double tmass = calc_xcm(top_x, index.size(), index.data(), top.atoms.atom, xcm, FALSE);
346 double Strans = calcTranslationalEntropy(tmass, T, P);
347 std::vector<gmx::RVec> x_com;
348 x_com.resize(top.atoms.nr);
349 for (int i = 0; i < top.atoms.nr; i++)
351 copy_rvec(top_x[i], x_com[i]);
353 (void)sub_xcm(as_rvec_array(x_com.data()), index.size(), index.data(), top.atoms.atom, xcm, FALSE);
357 principal_comp(index.size(), index.data(), top.atoms.atom, as_rvec_array(x_com.data()), trans, inertia);
358 bool linear = (inertia[XX] / inertia[YY] < linear_toler && inertia[XX] / inertia[ZZ] < linear_toler);
359 // (kJ/mol ps)^2/(Dalton nm^2 kJ/mol K) =
360 // KILO kg m^2 ps^2/(s^2 mol g/mol nm^2 K) =
361 // KILO^2 10^18 / 10^24 K = 1/K
362 double rot_const = gmx::square(PLANCK) / (8 * gmx::square(M_PI) * BOLTZ);
363 // Rotational temperature (1/K)
364 rvec theta = { 0, 0, 0 };
367 // For linear molecules the first element of the inertia
369 theta[0] = rot_const / inertia[1];
373 for (int m = 0; m < DIM; m++)
375 theta[m] = rot_const / inertia[m];
380 pr_rvec(debug, 0, "inertia", inertia, DIM, TRUE);
381 pr_rvec(debug, 0, "theta", theta, DIM, TRUE);
382 pr_rvecs(debug, 0, "trans", trans, DIM);
383 fprintf(debug, "linear molecule = %s\n", linear ? "true" : "no");
385 size_t nFreq = index.size() * DIM;
386 auto eFreq = gmx::arrayRefFromArray(eigfreq, nFreq);
387 double Svib = calcQuasiHarmonicEntropy(eFreq, T, linear, scale_factor);
389 double Srot = calcRotationalEntropy(T, top.atoms.nr, linear, theta, sigma_r);
390 fprintf(fp, "Translational entropy %g J/mol K\n", Strans);
391 fprintf(fp, "Rotational entropy %g J/mol K\n", Srot);
392 fprintf(fp, "Vibrational entropy %g J/mol K\n", Svib);
393 fprintf(fp, "Total entropy %g J/mol K\n", Svib + Strans + Srot);
394 double HeatCapacity = (calcTranslationalHeatCapacity() + calcRotationalHeatCapacity(linear)
395 + calcVibrationalHeatCapacity(eFreq, T, linear, scale_factor));
396 fprintf(fp, "Heat capacity %g J/mol K\n", HeatCapacity);
397 double Evib = (calcTranslationalInternalEnergy(T) + calcRotationalInternalEnergy(linear, T)
398 + calcVibrationalInternalEnergy(eFreq, T, linear, scale_factor));
399 fprintf(fp, "Internal energy %g kJ/mol\n", Evib);
400 double ZPE = calcZeroPointEnergy(eFreq, scale_factor);
401 fprintf(fp, "Zero-point energy %g kJ/mol\n", ZPE);
405 int gmx_nmeig(int argc, char* argv[])
407 const char* desc[] = {
408 "[THISMODULE] calculates the eigenvectors/values of a (Hessian) matrix,",
409 "which can be calculated with [gmx-mdrun].",
410 "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
411 "The structure is written first with t=0. The eigenvectors",
412 "are written as frames with the eigenvector number and eigenvalue",
413 "written as step number and timestamp, respectively.",
414 "The eigenvectors can be analyzed with [gmx-anaeig].",
415 "An ensemble of structures can be generated from the eigenvectors with",
416 "[gmx-nmens]. When mass weighting is used, the generated eigenvectors",
417 "will be scaled back to plain Cartesian coordinates before generating the",
418 "output. In this case, they will no longer be exactly orthogonal in the",
419 "standard Cartesian norm, but in the mass-weighted norm they would be.[PAR]",
420 "This program can be optionally used to compute quantum corrections to heat capacity",
421 "and enthalpy by providing an extra file argument [TT]-qcorr[tt]. See the GROMACS",
422 "manual, Chapter 1, for details. The result includes subtracting a harmonic",
423 "degree of freedom at the given temperature.",
424 "The total correction is printed on the terminal screen.",
425 "The recommended way of getting the corrections out is:[PAR]",
426 "[TT]gmx nmeig -s topol.tpr -f nm.mtx -first 7 -last 10000 -T 300 -qc [-constr][tt][PAR]",
427 "The [TT]-constr[tt] option should be used when bond constraints were used during the",
428 "simulation [BB]for all the covalent bonds[bb]. If this is not the case, ",
429 "you need to analyze the [TT]quant_corr.xvg[tt] file yourself.[PAR]",
430 "To make things more flexible, the program can also take virtual sites into account",
431 "when computing quantum corrections. When selecting [TT]-constr[tt] and",
432 "[TT]-qc[tt], the [TT]-begin[tt] and [TT]-end[tt] options will be set automatically as ",
434 "Based on a harmonic analysis of the normal mode frequencies,",
435 "thermochemical properties S0 (Standard Entropy),",
436 "Cv (Heat capacity at constant volume), Zero-point energy and the internal energy are",
437 "computed, much in the same manner as popular quantum chemistry",
441 static gmx_bool bM = TRUE, bCons = FALSE;
442 static int begin = 1, end = 50, maxspec = 4000, sigma_r = 1;
443 static real T = 298.15, width = 1, P = 1, scale_factor = 1;
444 static real linear_toler = 1e-5;
451 "Divide elements of Hessian by product of sqrt(mass) of involved "
452 "atoms prior to diagonalization. This should be used for 'Normal Modes' "
454 { "-first", FALSE, etINT, { &begin }, "First eigenvector to write away" },
459 "Last eigenvector to write away. -1 is use all dimensions." },
464 "Highest frequency (1/cm) to consider in the spectrum" },
469 "Temperature for computing entropy, quantum heat capacity and enthalpy when "
470 "using normal mode calculations to correct classical simulations" },
471 { "-P", FALSE, etREAL, { &P }, "Pressure (bar) when computing entropy" },
476 "Number of symmetric copies used when computing entropy. E.g. for water the "
477 "number is 2, for NH3 it is 3 and for methane it is 12." },
482 "Factor to scale frequencies before computing thermochemistry values" },
487 "Tolerance for determining whether a compound is linear as determined from the "
488 "ration of the moments inertion Ix/Iy and Ix/Iz." },
493 "If constraints were used in the simulation but not in the normal mode analysis "
494 "you will need to set this for computing the quantum corrections." },
499 "Width (sigma) of the gaussian peaks (1/cm) when generating a spectrum" }
501 FILE * out, *qc, *spec;
508 real qcvtot, qutot, qcv, qu;
510 real value, omega, nu;
511 real factor_gmx_to_omega2;
512 real factor_omega_to_wavenumber;
513 real* spectrum = nullptr;
515 gmx_output_env_t* oenv;
516 const char* qcleg[] = { "Heat Capacity cV (J/mol K)", "Enthalpy H (kJ/mol)" };
517 real* full_hessian = nullptr;
518 gmx_sparsematrix_t* sparse_hessian = nullptr;
521 { efMTX, "-f", "hessian", ffREAD }, { efTPR, nullptr, nullptr, ffREAD },
522 { efXVG, "-of", "eigenfreq", ffWRITE }, { efXVG, "-ol", "eigenval", ffWRITE },
523 { efXVG, "-os", "spectrum", ffOPTWR }, { efXVG, "-qc", "quant_corr", ffOPTWR },
524 { efTRN, "-v", "eigenvec", ffWRITE }
526 #define NFILE asize(fnm)
528 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
533 /* Read tpr file for volume and number of harmonic terms */
534 TpxFileHeader tpx = readTpxHeader(ftp2fn(efTPR, NFILE, fnm), true);
535 snew(top_x, tpx.natoms);
538 read_tpx(ftp2fn(efTPR, NFILE, fnm), nullptr, box, &natoms_tpx, top_x, nullptr, &mtop);
542 nharm = get_nharm(&mtop);
544 std::vector<int> atom_index = get_atom_index(&mtop);
546 top = gmx_mtop_t_to_t_topology(&mtop, true);
549 int ndim = DIM * atom_index.size();
551 if (opt2bSet("-qc", NFILE, fnm))
560 if (end == -1 || end > ndim)
564 printf("Using begin = %d and end = %d\n", begin, end);
566 /*open Hessian matrix */
568 gmx_mtxio_read(ftp2fn(efMTX, NFILE, fnm), &nrow, &ncol, &full_hessian, &sparse_hessian);
570 /* If the Hessian is in sparse format we can calculate max (ndim-1) eigenvectors,
571 * If this is not valid we convert to full matrix storage,
572 * but warn the user that we might run out of memory...
574 if ((sparse_hessian != nullptr) && (end == ndim))
576 fprintf(stderr, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
579 "Will try to allocate memory and convert to full matrix representation...\n");
581 size_t hessianSize = static_cast<size_t>(nrow) * static_cast<size_t>(ncol);
582 /* Allowing Hessians larger than INT_MAX probably only makes sense
583 * with (OpenMP) parallel diagonalization routines, since with a single
584 * thread it will takes months.
586 if (hessianSize > INT_MAX)
589 "Hessian size is %d x %d, which is larger than the maximum allowed %d "
591 nrow, ncol, INT_MAX);
593 snew(full_hessian, hessianSize);
594 for (i = 0; i < nrow * ncol; i++)
599 for (i = 0; i < sparse_hessian->nrow; i++)
601 for (j = 0; j < sparse_hessian->ndata[i]; j++)
603 k = sparse_hessian->data[i][j].col;
604 value = sparse_hessian->data[i][j].value;
605 full_hessian[i * ndim + k] = value;
606 full_hessian[k * ndim + i] = value;
609 gmx_sparsematrix_destroy(sparse_hessian);
610 sparse_hessian = nullptr;
611 fprintf(stderr, "Converted sparse to full matrix storage.\n");
614 snew(eigenvalues, nrow);
616 if (full_hessian != nullptr)
618 /* Using full matrix storage */
619 eigenvectors = allocateEigenvectors(nrow, begin, end, false);
621 nma_full_hessian(full_hessian, nrow, bM, &top, atom_index, begin, end, eigenvalues, eigenvectors);
625 assert(sparse_hessian);
626 /* Sparse memory storage, allocate memory for eigenvectors */
627 eigenvectors = allocateEigenvectors(nrow, begin, end, true);
629 nma_sparse_hessian(sparse_hessian, bM, &top, atom_index, end, eigenvalues, eigenvectors);
632 /* check the output, first 6 eigenvalues should be reasonably small */
633 gmx_bool bSuck = FALSE;
634 for (i = begin - 1; (i < 6); i++)
636 if (std::abs(eigenvalues[i]) > 1.0e-3)
643 fprintf(stderr, "\nOne of the lowest 6 eigenvalues has a non-zero value.\n");
644 fprintf(stderr, "This could mean that the reference structure was not\n");
645 fprintf(stderr, "properly energy minimized.\n");
648 /* now write the output */
649 fprintf(stderr, "Writing eigenvalues...\n");
650 out = xvgropen(opt2fn("-ol", NFILE, fnm), "Eigenvalues", "Eigenvalue index",
651 "Eigenvalue [Gromacs units]", oenv);
652 if (output_env_get_print_xvgr_codes(oenv))
656 fprintf(out, "@ subtitle \"mass weighted\"\n");
660 fprintf(out, "@ subtitle \"not mass weighted\"\n");
664 for (i = 0; i <= (end - begin); i++)
666 fprintf(out, "%6d %15g\n", begin + i, eigenvalues[i]);
671 if (opt2bSet("-qc", NFILE, fnm))
673 qc = xvgropen(opt2fn("-qc", NFILE, fnm), "Quantum Corrections", "Eigenvector index", "", oenv);
674 xvgr_legend(qc, asize(qcleg), qcleg, oenv);
681 printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
683 out = xvgropen(opt2fn("-of", NFILE, fnm), "Eigenfrequencies", "Eigenvector index",
684 "Wavenumber [cm\\S-1\\N]", oenv);
685 if (output_env_get_print_xvgr_codes(oenv))
689 fprintf(out, "@ subtitle \"mass weighted\"\n");
693 fprintf(out, "@ subtitle \"not mass weighted\"\n");
698 if (opt2bSet("-os", NFILE, fnm) && (maxspec > 0))
700 snew(spectrum, maxspec);
701 spec = xvgropen(opt2fn("-os", NFILE, fnm),
702 "Vibrational spectrum based on harmonic approximation",
703 "\\f{12}w\\f{4} (cm\\S-1\\N)", "Intensity [Gromacs units]", oenv);
704 for (i = 0; (i < maxspec); i++)
710 /* Gromacs units are kJ/(mol*nm*nm*amu),
711 * where amu is the atomic mass unit.
713 * For the eigenfrequencies we want to convert this to spectroscopic absorption
714 * wavenumbers given in cm^(-1), which is the frequency divided by the speed of
715 * light. Do this by first converting to omega^2 (units 1/s), take the square
716 * root, and finally divide by the speed of light (nm/ps in gromacs).
718 factor_gmx_to_omega2 = 1.0E21 / (AVOGADRO * AMU);
719 factor_omega_to_wavenumber = 1.0E-5 / (2.0 * M_PI * SPEED_OF_LIGHT);
722 for (i = begin; (i <= end); i++)
724 value = eigenvalues[i - begin];
729 omega = std::sqrt(value * factor_gmx_to_omega2);
730 nu = 1e-12 * omega / (2 * M_PI);
731 value = omega * factor_omega_to_wavenumber;
732 fprintf(out, "%6d %15g\n", i, value);
735 wfac = eigenvalues[i - begin] / (width * std::sqrt(2 * M_PI));
736 for (j = 0; (j < maxspec); j++)
738 spectrum[j] += wfac * std::exp(-gmx::square(j - value) / (2 * gmx::square(width)));
743 qcv = cv_corr(nu, T);
750 fprintf(qc, "%6d %15g %15g\n", i, qcv, qu);
756 if (value >= maxspec)
758 printf("WARNING: high frequencies encountered (%g cm^-1).\n", value);
759 printf("Your calculations may be incorrect due to e.g. improper minimization of\n");
760 printf("your starting structure or due to issues in your topology.\n");
764 for (j = 0; (j < maxspec); j++)
766 fprintf(spec, "%10g %10g\n", 1.0 * j, spectrum[j]);
772 printf("Quantum corrections for harmonic degrees of freedom\n");
773 printf("Use appropriate -first and -last options to get reliable results.\n");
774 printf("There were %d constraints in the simulation\n", nharm);
775 printf("Total correction to cV = %g J/mol K\n", qcvtot);
776 printf("Total correction to H = %g kJ/mol\n", qutot);
778 please_cite(stdout, "Caleman2011b");
780 /* Writing eigenvectors. Note that if mass scaling was used, the eigenvectors
781 * were scaled back from mass weighted cartesian to plain cartesian in the
782 * nma_full_hessian() or nma_sparse_hessian() routines. Mass scaled vectors
783 * will not be strictly orthogonal in plain cartesian scalar products.
785 const real* eigenvectorPtr;
786 if (full_hessian != nullptr)
788 eigenvectorPtr = eigenvectors;
792 /* The sparse matrix diagonalization store all eigenvectors up to end */
793 eigenvectorPtr = eigenvectors + (begin - 1) * atom_index.size();
795 write_eigenvectors(opt2fn("-v", NFILE, fnm), atom_index.size(), eigenvectorPtr, FALSE, begin,
796 end, eWXR_NO, nullptr, FALSE, top_x, bM, eigenvalues);
800 analyzeThermochemistry(stdout, top, top_x, atom_index, eigenvalues, T, P, sigma_r,
801 scale_factor, linear_toler);
802 please_cite(stdout, "Spoel2018a");
806 printf("Cannot compute entropy when -first = %d\n", begin);