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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/mtxio.h"
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/eigio.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/gmxana/gstat.h"
53 #include "gromacs/gmxana/princ.h"
54 #include "gromacs/linearalgebra/eigensolver.h"
55 #include "gromacs/linearalgebra/sparsematrix.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vecdump.h"
60 #include "gromacs/topology/ifunc.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/gmxassert.h"
67 #include "gromacs/utility/pleasecite.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "thermochemistry.h"
72 static double cv_corr(double nu
, double T
)
74 double x
= PLANCK
* nu
/ (BOLTZ
* T
);
75 double ex
= std::exp(x
);
83 return BOLTZ
* KILO
* (ex
* gmx::square(x
) / gmx::square(ex
- 1) - 1);
87 static double u_corr(double nu
, double T
)
89 double x
= PLANCK
* nu
/ (BOLTZ
* T
);
90 double ex
= std::exp(x
);
98 return BOLTZ
* T
* (0.5 * x
- 1 + x
/ (ex
- 1));
102 static size_t get_nharm_mt(const gmx_moltype_t
* mt
)
104 static int harm_func
[] = { F_BONDS
};
108 for (i
= 0; (i
< asize(harm_func
)); i
++)
111 nh
+= mt
->ilist
[ft
].size() / (interaction_function
[ft
].nratoms
+ 1);
116 static int get_nharm(const gmx_mtop_t
* mtop
)
120 for (const gmx_molblock_t
& molb
: mtop
->molblock
)
122 nh
+= molb
.nmol
* get_nharm_mt(&(mtop
->moltype
[molb
.type
]));
127 static void nma_full_hessian(real
* hess
,
130 const t_topology
* top
,
131 gmx::ArrayRef
<const int> atom_index
,
139 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
143 for (int i
= 0; (i
< atom_index
.ssize()); i
++)
145 size_t ai
= atom_index
[i
];
146 for (size_t j
= 0; (j
< DIM
); j
++)
148 for (int k
= 0; (k
< atom_index
.ssize()); k
++)
150 size_t ak
= atom_index
[k
];
151 mass_fac
= gmx::invsqrt(top
->atoms
.atom
[ai
].m
* top
->atoms
.atom
[ak
].m
);
152 for (size_t l
= 0; (l
< DIM
); l
++)
154 hess
[(i
* DIM
+ j
) * ndim
+ k
* DIM
+ l
] *= mass_fac
;
161 /* call diagonalization routine. */
163 fprintf(stderr
, "\nDiagonalizing to find vectors %d through %d...\n", begin
, end
);
166 eigensolver(hess
, ndim
, begin
- 1, end
- 1, eigenvalues
, eigenvectors
);
168 /* And scale the output eigenvectors */
169 if (bM
&& eigenvectors
!= nullptr)
171 for (int i
= 0; i
< (end
- begin
+ 1); i
++)
173 for (gmx::index j
= 0; j
< atom_index
.ssize(); j
++)
175 size_t aj
= atom_index
[j
];
176 mass_fac
= gmx::invsqrt(top
->atoms
.atom
[aj
].m
);
177 for (size_t k
= 0; (k
< DIM
); k
++)
179 eigenvectors
[i
* ndim
+ j
* DIM
+ k
] *= mass_fac
;
187 static void nma_sparse_hessian(gmx_sparsematrix_t
* sparse_hessian
,
189 const t_topology
* top
,
190 gmx::ArrayRef
<const int> 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,
207 "NULL matrix pointer provided to nma_sparse_hessian");
211 for (int iatom
= 0; (iatom
< atom_index
.ssize()); iatom
++)
213 size_t ai
= atom_index
[iatom
];
214 for (size_t j
= 0; (j
< DIM
); j
++)
216 row
= DIM
* iatom
+ j
;
217 for (k
= 0; k
< sparse_hessian
->ndata
[row
]; k
++)
219 col
= sparse_hessian
->data
[row
][k
].col
;
221 size_t ak
= atom_index
[katom
];
222 mass_fac
= gmx::invsqrt(top
->atoms
.atom
[ai
].m
* top
->atoms
.atom
[ak
].m
);
223 sparse_hessian
->data
[row
][k
].value
*= mass_fac
;
228 fprintf(stderr
, "\nDiagonalizing to find eigenvectors 1 through %d...\n", neig
);
231 sparse_eigensolver(sparse_hessian
, neig
, eigenvalues
, eigenvectors
, 10000000);
233 /* Scale output eigenvectors */
234 if (bM
&& eigenvectors
!= nullptr)
236 for (i
= 0; i
< neig
; i
++)
238 for (gmx::index j
= 0; j
< atom_index
.ssize(); j
++)
240 size_t aj
= atom_index
[j
];
241 mass_fac
= gmx::invsqrt(top
->atoms
.atom
[aj
].m
);
242 for (k
= 0; (k
< DIM
); k
++)
244 eigenvectors
[i
* ndim
+ j
* DIM
+ k
] *= mass_fac
;
252 /* Returns a pointer for eigenvector storage */
253 static real
* allocateEigenvectors(int nrow
, int first
, int last
, bool ignoreBegin
)
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
)
272 "You asked to store %d eigenvectors of size %d, which requires more than the "
273 "supported %d elements; %sdecrease -last",
274 numVector
, nrow
, INT_MAX
, ignoreBegin
? "" : "increase -first and/or ");
278 snew(eigenvectors
, vectorsSize
);
283 /*! \brief Compute heat capacity due to translational motion
285 static double calcTranslationalHeatCapacity()
290 /*! \brief Compute internal energy due to translational motion
292 static double calcTranslationalInternalEnergy(double T
)
294 return BOLTZ
* T
* 1.5;
297 /*! \brief Compute heat capacity due to rotational motion
299 * \param[in] linear Should be TRUE if this is a linear molecule
300 * \param[in] T Temperature
301 * \return The heat capacity at constant volume
303 static double calcRotationalInternalEnergy(gmx_bool linear
, double T
)
311 return BOLTZ
* T
* 1.5;
315 /*! \brief Compute heat capacity due to rotational motion
317 * \param[in] linear Should be TRUE if this is a linear molecule
318 * \return The heat capacity at constant volume
320 static double calcRotationalHeatCapacity(gmx_bool linear
)
332 static void analyzeThermochemistry(FILE* fp
,
333 const t_topology
& top
,
335 gmx::ArrayRef
<const int> atom_index
,
343 std::vector
<int> index(atom_index
.begin(), atom_index
.end());
346 double tmass
= calc_xcm(top_x
, index
.size(), index
.data(), top
.atoms
.atom
, xcm
, FALSE
);
347 double Strans
= calcTranslationalEntropy(tmass
, T
, P
);
348 std::vector
<gmx::RVec
> x_com
;
349 x_com
.resize(top
.atoms
.nr
);
350 for (int i
= 0; i
< top
.atoms
.nr
; i
++)
352 copy_rvec(top_x
[i
], x_com
[i
]);
354 (void)sub_xcm(as_rvec_array(x_com
.data()), index
.size(), index
.data(), top
.atoms
.atom
, xcm
, FALSE
);
358 principal_comp(index
.size(), index
.data(), top
.atoms
.atom
, as_rvec_array(x_com
.data()), trans
, inertia
);
359 bool linear
= (inertia
[XX
] / inertia
[YY
] < linear_toler
&& inertia
[XX
] / inertia
[ZZ
] < linear_toler
);
360 // (kJ/mol ps)^2/(Dalton nm^2 kJ/mol K) =
361 // KILO kg m^2 ps^2/(s^2 mol g/mol nm^2 K) =
362 // KILO^2 10^18 / 10^24 K = 1/K
363 double rot_const
= gmx::square(PLANCK
) / (8 * gmx::square(M_PI
) * BOLTZ
);
364 // Rotational temperature (1/K)
365 rvec theta
= { 0, 0, 0 };
368 // For linear molecules the first element of the inertia
370 theta
[0] = rot_const
/ inertia
[1];
374 for (int m
= 0; m
< DIM
; m
++)
376 theta
[m
] = rot_const
/ inertia
[m
];
381 pr_rvec(debug
, 0, "inertia", inertia
, DIM
, TRUE
);
382 pr_rvec(debug
, 0, "theta", theta
, DIM
, TRUE
);
383 pr_rvecs(debug
, 0, "trans", trans
, DIM
);
384 fprintf(debug
, "linear molecule = %s\n", linear
? "true" : "no");
386 size_t nFreq
= index
.size() * DIM
;
387 auto eFreq
= gmx::arrayRefFromArray(eigfreq
, nFreq
);
388 double Svib
= calcQuasiHarmonicEntropy(eFreq
, T
, linear
, scale_factor
);
390 double Srot
= calcRotationalEntropy(T
, top
.atoms
.nr
, linear
, theta
, sigma_r
);
391 fprintf(fp
, "Translational entropy %g J/mol K\n", Strans
);
392 fprintf(fp
, "Rotational entropy %g J/mol K\n", Srot
);
393 fprintf(fp
, "Vibrational entropy %g J/mol K\n", Svib
);
394 fprintf(fp
, "Total entropy %g J/mol K\n", Svib
+ Strans
+ Srot
);
395 double HeatCapacity
= (calcTranslationalHeatCapacity() + calcRotationalHeatCapacity(linear
)
396 + calcVibrationalHeatCapacity(eFreq
, T
, linear
, scale_factor
));
397 fprintf(fp
, "Heat capacity %g J/mol K\n", HeatCapacity
);
398 double Evib
= (calcTranslationalInternalEnergy(T
) + calcRotationalInternalEnergy(linear
, T
)
399 + calcVibrationalInternalEnergy(eFreq
, T
, linear
, scale_factor
));
400 fprintf(fp
, "Internal energy %g kJ/mol\n", Evib
);
401 double ZPE
= calcZeroPointEnergy(eFreq
, scale_factor
);
402 fprintf(fp
, "Zero-point energy %g kJ/mol\n", ZPE
);
406 int gmx_nmeig(int argc
, char* argv
[])
408 const char* desc
[] = {
409 "[THISMODULE] calculates the eigenvectors/values of a (Hessian) matrix,",
410 "which can be calculated with [gmx-mdrun].",
411 "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
412 "The structure is written first with t=0. The eigenvectors",
413 "are written as frames with the eigenvector number and eigenvalue",
414 "written as step number and timestamp, respectively.",
415 "The eigenvectors can be analyzed with [gmx-anaeig].",
416 "An ensemble of structures can be generated from the eigenvectors with",
417 "[gmx-nmens]. When mass weighting is used, the generated eigenvectors",
418 "will be scaled back to plain Cartesian coordinates before generating the",
419 "output. In this case, they will no longer be exactly orthogonal in the",
420 "standard Cartesian norm, but in the mass-weighted norm they would be.[PAR]",
421 "This program can be optionally used to compute quantum corrections to heat capacity",
422 "and enthalpy by providing an extra file argument [TT]-qcorr[tt]. See the GROMACS",
423 "manual, Chapter 1, for details. The result includes subtracting a harmonic",
424 "degree of freedom at the given temperature.",
425 "The total correction is printed on the terminal screen.",
426 "The recommended way of getting the corrections out is:[PAR]",
427 "[TT]gmx nmeig -s topol.tpr -f nm.mtx -first 7 -last 10000 -T 300 -qc [-constr][tt][PAR]",
428 "The [TT]-constr[tt] option should be used when bond constraints were used during the",
429 "simulation [BB]for all the covalent bonds[bb]. If this is not the case, ",
430 "you need to analyze the [TT]quant_corr.xvg[tt] file yourself.[PAR]",
431 "To make things more flexible, the program can also take virtual sites into account",
432 "when computing quantum corrections. When selecting [TT]-constr[tt] and",
433 "[TT]-qc[tt], the [TT]-begin[tt] and [TT]-end[tt] options will be set automatically as ",
435 "Based on a harmonic analysis of the normal mode frequencies,",
436 "thermochemical properties S0 (Standard Entropy),",
437 "Cv (Heat capacity at constant volume), Zero-point energy and the internal energy are",
438 "computed, much in the same manner as popular quantum chemistry",
442 static gmx_bool bM
= TRUE
, bCons
= FALSE
;
443 static int begin
= 1, end
= 50, maxspec
= 4000, sigma_r
= 1;
444 static real T
= 298.15, width
= 1, P
= 1, scale_factor
= 1;
445 static real linear_toler
= 1e-5;
452 "Divide elements of Hessian by product of sqrt(mass) of involved "
453 "atoms prior to diagonalization. This should be used for 'Normal Modes' "
455 { "-first", FALSE
, etINT
, { &begin
}, "First eigenvector to write away" },
460 "Last eigenvector to write away. -1 is use all dimensions." },
465 "Highest frequency (1/cm) to consider in the spectrum" },
470 "Temperature for computing entropy, quantum heat capacity and enthalpy when "
471 "using normal mode calculations to correct classical simulations" },
472 { "-P", FALSE
, etREAL
, { &P
}, "Pressure (bar) when computing entropy" },
477 "Number of symmetric copies used when computing entropy. E.g. for water the "
478 "number is 2, for NH3 it is 3 and for methane it is 12." },
483 "Factor to scale frequencies before computing thermochemistry values" },
488 "Tolerance for determining whether a compound is linear as determined from the "
489 "ration of the moments inertion Ix/Iy and Ix/Iz." },
494 "If constraints were used in the simulation but not in the normal mode analysis "
495 "you will need to set this for computing the quantum corrections." },
500 "Width (sigma) of the gaussian peaks (1/cm) when generating a spectrum" }
502 FILE * out
, *qc
, *spec
;
509 real qcvtot
, qutot
, qcv
, qu
;
511 real value
, omega
, nu
;
512 real factor_gmx_to_omega2
;
513 real factor_omega_to_wavenumber
;
514 real
* spectrum
= nullptr;
516 gmx_output_env_t
* oenv
;
517 const char* qcleg
[] = { "Heat Capacity cV (J/mol K)", "Enthalpy H (kJ/mol)" };
518 real
* full_hessian
= nullptr;
519 gmx_sparsematrix_t
* sparse_hessian
= nullptr;
522 { efMTX
, "-f", "hessian", ffREAD
}, { efTPR
, nullptr, nullptr, ffREAD
},
523 { efXVG
, "-of", "eigenfreq", ffWRITE
}, { efXVG
, "-ol", "eigenval", ffWRITE
},
524 { efXVG
, "-os", "spectrum", ffOPTWR
}, { efXVG
, "-qc", "quant_corr", ffOPTWR
},
525 { efTRN
, "-v", "eigenvec", ffWRITE
}
527 #define NFILE asize(fnm)
529 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
, asize(desc
), desc
, 0, nullptr, &oenv
))
534 /* Read tpr file for volume and number of harmonic terms */
535 TpxFileHeader tpx
= readTpxHeader(ftp2fn(efTPR
, NFILE
, fnm
), true);
536 snew(top_x
, tpx
.natoms
);
539 read_tpx(ftp2fn(efTPR
, NFILE
, fnm
), nullptr, box
, &natoms_tpx
, top_x
, nullptr, &mtop
);
543 nharm
= get_nharm(&mtop
);
545 std::vector
<int> atom_index
= get_atom_index(&mtop
);
547 top
= gmx_mtop_t_to_t_topology(&mtop
, true);
550 int ndim
= DIM
* atom_index
.size();
552 if (opt2bSet("-qc", NFILE
, fnm
))
561 if (end
== -1 || end
> ndim
)
565 printf("Using begin = %d and end = %d\n", begin
, end
);
567 /*open Hessian matrix */
569 gmx_mtxio_read(ftp2fn(efMTX
, NFILE
, fnm
), &nrow
, &ncol
, &full_hessian
, &sparse_hessian
);
571 /* If the Hessian is in sparse format we can calculate max (ndim-1) eigenvectors,
572 * If this is not valid we convert to full matrix storage,
573 * but warn the user that we might run out of memory...
575 if ((sparse_hessian
!= nullptr) && (end
== ndim
))
577 fprintf(stderr
, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
580 "Will try to allocate memory and convert to full matrix representation...\n");
582 size_t hessianSize
= static_cast<size_t>(nrow
) * static_cast<size_t>(ncol
);
583 /* Allowing Hessians larger than INT_MAX probably only makes sense
584 * with (OpenMP) parallel diagonalization routines, since with a single
585 * thread it will takes months.
587 if (hessianSize
> INT_MAX
)
590 "Hessian size is %d x %d, which is larger than the maximum allowed %d "
592 nrow
, ncol
, INT_MAX
);
594 snew(full_hessian
, hessianSize
);
595 for (i
= 0; i
< nrow
* ncol
; i
++)
600 for (i
= 0; i
< sparse_hessian
->nrow
; i
++)
602 for (j
= 0; j
< sparse_hessian
->ndata
[i
]; j
++)
604 k
= sparse_hessian
->data
[i
][j
].col
;
605 value
= sparse_hessian
->data
[i
][j
].value
;
606 full_hessian
[i
* ndim
+ k
] = value
;
607 full_hessian
[k
* ndim
+ i
] = value
;
610 gmx_sparsematrix_destroy(sparse_hessian
);
611 sparse_hessian
= nullptr;
612 fprintf(stderr
, "Converted sparse to full matrix storage.\n");
615 snew(eigenvalues
, nrow
);
617 if (full_hessian
!= nullptr)
619 /* Using full matrix storage */
620 eigenvectors
= allocateEigenvectors(nrow
, begin
, end
, false);
622 nma_full_hessian(full_hessian
, nrow
, bM
, &top
, atom_index
, begin
, end
, eigenvalues
, eigenvectors
);
626 assert(sparse_hessian
);
627 /* Sparse memory storage, allocate memory for eigenvectors */
628 eigenvectors
= allocateEigenvectors(nrow
, begin
, end
, true);
630 nma_sparse_hessian(sparse_hessian
, bM
, &top
, atom_index
, end
, eigenvalues
, eigenvectors
);
633 /* check the output, first 6 eigenvalues should be reasonably small */
634 gmx_bool bSuck
= FALSE
;
635 for (i
= begin
- 1; (i
< 6); i
++)
637 if (std::abs(eigenvalues
[i
]) > 1.0e-3)
644 fprintf(stderr
, "\nOne of the lowest 6 eigenvalues has a non-zero value.\n");
645 fprintf(stderr
, "This could mean that the reference structure was not\n");
646 fprintf(stderr
, "properly energy minimized.\n");
649 /* now write the output */
650 fprintf(stderr
, "Writing eigenvalues...\n");
651 out
= xvgropen(opt2fn("-ol", NFILE
, fnm
), "Eigenvalues", "Eigenvalue index",
652 "Eigenvalue [Gromacs units]", oenv
);
653 if (output_env_get_print_xvgr_codes(oenv
))
657 fprintf(out
, "@ subtitle \"mass weighted\"\n");
661 fprintf(out
, "@ subtitle \"not mass weighted\"\n");
665 for (i
= 0; i
<= (end
- begin
); i
++)
667 fprintf(out
, "%6d %15g\n", begin
+ i
, eigenvalues
[i
]);
672 if (opt2bSet("-qc", NFILE
, fnm
))
674 qc
= xvgropen(opt2fn("-qc", NFILE
, fnm
), "Quantum Corrections", "Eigenvector index", "", oenv
);
675 xvgr_legend(qc
, asize(qcleg
), qcleg
, oenv
);
682 printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
684 out
= xvgropen(opt2fn("-of", NFILE
, fnm
), "Eigenfrequencies", "Eigenvector index",
685 "Wavenumber [cm\\S-1\\N]", oenv
);
686 if (output_env_get_print_xvgr_codes(oenv
))
690 fprintf(out
, "@ subtitle \"mass weighted\"\n");
694 fprintf(out
, "@ subtitle \"not mass weighted\"\n");
699 if (opt2bSet("-os", NFILE
, fnm
) && (maxspec
> 0))
701 snew(spectrum
, maxspec
);
702 spec
= xvgropen(opt2fn("-os", NFILE
, fnm
),
703 "Vibrational spectrum based on harmonic approximation",
704 "\\f{12}w\\f{4} (cm\\S-1\\N)", "Intensity [Gromacs units]", oenv
);
705 for (i
= 0; (i
< maxspec
); i
++)
711 /* Gromacs units are kJ/(mol*nm*nm*amu),
712 * where amu is the atomic mass unit.
714 * For the eigenfrequencies we want to convert this to spectroscopic absorption
715 * wavenumbers given in cm^(-1), which is the frequency divided by the speed of
716 * light. Do this by first converting to omega^2 (units 1/s), take the square
717 * root, and finally divide by the speed of light (nm/ps in gromacs).
719 factor_gmx_to_omega2
= 1.0E21
/ (AVOGADRO
* AMU
);
720 factor_omega_to_wavenumber
= 1.0E-5 / (2.0 * M_PI
* SPEED_OF_LIGHT
);
723 for (i
= begin
; (i
<= end
); i
++)
725 value
= eigenvalues
[i
- begin
];
730 omega
= std::sqrt(value
* factor_gmx_to_omega2
);
731 nu
= 1e-12 * omega
/ (2 * M_PI
);
732 value
= omega
* factor_omega_to_wavenumber
;
733 fprintf(out
, "%6d %15g\n", i
, value
);
736 wfac
= eigenvalues
[i
- begin
] / (width
* std::sqrt(2 * M_PI
));
737 for (j
= 0; (j
< maxspec
); j
++)
739 spectrum
[j
] += wfac
* std::exp(-gmx::square(j
- value
) / (2 * gmx::square(width
)));
744 qcv
= cv_corr(nu
, T
);
751 fprintf(qc
, "%6d %15g %15g\n", i
, qcv
, qu
);
757 if (value
>= maxspec
)
759 printf("WARNING: high frequencies encountered (%g cm^-1).\n", value
);
760 printf("Your calculations may be incorrect due to e.g. improper minimization of\n");
761 printf("your starting structure or due to issues in your topology.\n");
765 for (j
= 0; (j
< maxspec
); j
++)
767 fprintf(spec
, "%10g %10g\n", 1.0 * j
, spectrum
[j
]);
773 printf("Quantum corrections for harmonic degrees of freedom\n");
774 printf("Use appropriate -first and -last options to get reliable results.\n");
775 printf("There were %d constraints in the simulation\n", nharm
);
776 printf("Total correction to cV = %g J/mol K\n", qcvtot
);
777 printf("Total correction to H = %g kJ/mol\n", qutot
);
779 please_cite(stdout
, "Caleman2011b");
781 /* Writing eigenvectors. Note that if mass scaling was used, the eigenvectors
782 * were scaled back from mass weighted cartesian to plain cartesian in the
783 * nma_full_hessian() or nma_sparse_hessian() routines. Mass scaled vectors
784 * will not be strictly orthogonal in plain cartesian scalar products.
786 const real
* eigenvectorPtr
;
787 if (full_hessian
!= nullptr)
789 eigenvectorPtr
= eigenvectors
;
793 /* The sparse matrix diagonalization store all eigenvectors up to end */
794 eigenvectorPtr
= eigenvectors
+ (begin
- 1) * atom_index
.size();
796 write_eigenvectors(opt2fn("-v", NFILE
, fnm
), atom_index
.size(), eigenvectorPtr
, FALSE
, begin
,
797 end
, eWXR_NO
, nullptr, FALSE
, top_x
, bM
, eigenvalues
);
801 analyzeThermochemistry(stdout
, top
, top_x
, atom_index
, eigenvalues
, T
, P
, sigma_r
,
802 scale_factor
, linear_toler
);
803 please_cite(stdout
, "Spoel2018a");
807 printf("Cannot compute entropy when -first = %d\n", begin
);