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, 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.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/copyrite.h"
44 #include "gromacs/fileio/mtxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/eigio.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/legacyheaders/txtdump.h"
50 #include "gromacs/legacyheaders/types/ifunc.h"
51 #include "gromacs/linearalgebra/eigensolver.h"
52 #include "gromacs/linearalgebra/sparsematrix.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/smalloc.h"
62 static double cv_corr(double nu
, double T
)
64 double x
= PLANCK
*nu
/(BOLTZ
*T
);
65 double ex
= std::exp(x
);
73 return BOLTZ
*KILO
*(ex
*sqr(x
)/sqr(ex
-1) - 1);
77 static double u_corr(double nu
, double T
)
79 double x
= PLANCK
*nu
/(BOLTZ
*T
);
80 double ex
= std::exp(x
);
88 return BOLTZ
*T
*(0.5*x
- 1 + x
/(ex
-1));
92 static int get_nharm_mt(gmx_moltype_t
*mt
)
94 static int harm_func
[] = { F_BONDS
};
98 for (i
= 0; (i
< asize(harm_func
)); i
++)
101 nh
+= mt
->ilist
[ft
].nr
/(interaction_function
[ft
].nratoms
+1);
106 static int get_nvsite_mt(gmx_moltype_t
*mt
)
108 static int vs_func
[] = {
109 F_VSITE2
, F_VSITE3
, F_VSITE3FD
, F_VSITE3FAD
,
110 F_VSITE3OUT
, F_VSITE4FD
, F_VSITE4FDN
, F_VSITEN
115 for (i
= 0; (i
< asize(vs_func
)); i
++)
118 nh
+= mt
->ilist
[ft
].nr
/(interaction_function
[ft
].nratoms
+1);
123 static int get_nharm(gmx_mtop_t
*mtop
, int *nvsites
)
129 for (j
= 0; (j
< mtop
->nmolblock
); j
++)
131 mt
= mtop
->molblock
[j
].type
;
132 nh
+= mtop
->molblock
[j
].nmol
* get_nharm_mt(&(mtop
->moltype
[mt
]));
133 nv
+= mtop
->molblock
[j
].nmol
* get_nvsite_mt(&(mtop
->moltype
[mt
]));
140 nma_full_hessian(real
* hess
,
153 natoms
= top
->atoms
.nr
;
155 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
159 for (i
= 0; (i
< natoms
); i
++)
161 for (j
= 0; (j
< DIM
); j
++)
163 for (k
= 0; (k
< natoms
); k
++)
165 mass_fac
= gmx_invsqrt(top
->atoms
.atom
[i
].m
*top
->atoms
.atom
[k
].m
);
166 for (l
= 0; (l
< DIM
); l
++)
168 hess
[(i
*DIM
+j
)*ndim
+k
*DIM
+l
] *= mass_fac
;
175 /* call diagonalization routine. */
177 fprintf(stderr
, "\nDiagonalizing to find vectors %d through %d...\n", begin
, end
);
180 eigensolver(hess
, ndim
, begin
-1, end
-1, eigenvalues
, eigenvectors
);
182 /* And scale the output eigenvectors */
183 if (bM
&& eigenvectors
!= NULL
)
185 for (i
= 0; i
< (end
-begin
+1); i
++)
187 for (j
= 0; j
< natoms
; j
++)
189 mass_fac
= gmx_invsqrt(top
->atoms
.atom
[j
].m
);
190 for (k
= 0; (k
< DIM
); k
++)
192 eigenvectors
[i
*ndim
+j
*DIM
+k
] *= mass_fac
;
202 nma_sparse_hessian(gmx_sparsematrix_t
* sparse_hessian
,
216 natoms
= top
->atoms
.nr
;
219 /* Cannot check symmetry since we only store half matrix */
220 /* divide elements hess[i][j] by sqrt(mas[i])*sqrt(mas[j]) when required */
222 GMX_RELEASE_ASSERT(sparse_hessian
!= NULL
, "NULL matrix pointer provided to nma_sparse_hessian");
226 for (iatom
= 0; (iatom
< natoms
); iatom
++)
228 for (j
= 0; (j
< DIM
); j
++)
231 for (k
= 0; k
< sparse_hessian
->ndata
[row
]; k
++)
233 col
= sparse_hessian
->data
[row
][k
].col
;
235 mass_fac
= gmx_invsqrt(top
->atoms
.atom
[iatom
].m
*top
->atoms
.atom
[katom
].m
);
236 sparse_hessian
->data
[row
][k
].value
*= mass_fac
;
241 fprintf(stderr
, "\nDiagonalizing to find eigenvectors 1 through %d...\n", neig
);
244 sparse_eigensolver(sparse_hessian
, neig
, eigenvalues
, eigenvectors
, 10000000);
246 /* Scale output eigenvectors */
247 if (bM
&& eigenvectors
!= NULL
)
249 for (i
= 0; i
< neig
; i
++)
251 for (j
= 0; j
< natoms
; j
++)
253 mass_fac
= gmx_invsqrt(top
->atoms
.atom
[j
].m
);
254 for (k
= 0; (k
< DIM
); k
++)
256 eigenvectors
[i
*ndim
+j
*DIM
+k
] *= mass_fac
;
265 int gmx_nmeig(int argc
, char *argv
[])
267 const char *desc
[] = {
268 "[THISMODULE] calculates the eigenvectors/values of a (Hessian) matrix,",
269 "which can be calculated with [gmx-mdrun].",
270 "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
271 "The structure is written first with t=0. The eigenvectors",
272 "are written as frames with the eigenvector number as timestamp.",
273 "The eigenvectors can be analyzed with [gmx-anaeig].",
274 "An ensemble of structures can be generated from the eigenvectors with",
275 "[gmx-nmens]. When mass weighting is used, the generated eigenvectors",
276 "will be scaled back to plain Cartesian coordinates before generating the",
277 "output. In this case, they will no longer be exactly orthogonal in the",
278 "standard Cartesian norm, but in the mass-weighted norm they would be.[PAR]",
279 "This program can be optionally used to compute quantum corrections to heat capacity",
280 "and enthalpy by providing an extra file argument [TT]-qcorr[tt]. See the GROMACS",
281 "manual, Chapter 1, for details. The result includes subtracting a harmonic",
282 "degree of freedom at the given temperature.",
283 "The total correction is printed on the terminal screen.",
284 "The recommended way of getting the corrections out is:[PAR]",
285 "[TT]gmx nmeig -s topol.tpr -f nm.mtx -first 7 -last 10000 -T 300 -qc [-constr][tt][PAR]",
286 "The [TT]-constr[tt] option should be used when bond constraints were used during the",
287 "simulation [BB]for all the covalent bonds[bb]. If this is not the case, ",
288 "you need to analyze the [TT]quant_corr.xvg[tt] file yourself.[PAR]",
289 "To make things more flexible, the program can also take virtual sites into account",
290 "when computing quantum corrections. When selecting [TT]-constr[tt] and",
291 "[TT]-qc[tt], the [TT]-begin[tt] and [TT]-end[tt] options will be set automatically as well.",
292 "Again, if you think you know it better, please check the [TT]eigenfreq.xvg[tt]",
296 static gmx_bool bM
= TRUE
, bCons
= FALSE
;
297 static int begin
= 1, end
= 50, maxspec
= 4000;
298 static real T
= 298.15, width
= 1;
301 { "-m", FALSE
, etBOOL
, {&bM
},
302 "Divide elements of Hessian by product of sqrt(mass) of involved "
303 "atoms prior to diagonalization. This should be used for 'Normal Modes' "
305 { "-first", FALSE
, etINT
, {&begin
},
306 "First eigenvector to write away" },
307 { "-last", FALSE
, etINT
, {&end
},
308 "Last eigenvector to write away" },
309 { "-maxspec", FALSE
, etINT
, {&maxspec
},
310 "Highest frequency (1/cm) to consider in the spectrum" },
311 { "-T", FALSE
, etREAL
, {&T
},
312 "Temperature for computing quantum heat capacity and enthalpy when using normal mode calculations to correct classical simulations" },
313 { "-constr", FALSE
, etBOOL
, {&bCons
},
314 "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." },
315 { "-width", FALSE
, etREAL
, {&width
},
316 "Width (sigma) of the gaussian peaks (1/cm) when generating a spectrum" }
318 FILE *out
, *qc
, *spec
;
325 real qcvtot
, qutot
, qcv
, qu
;
326 int natoms
, ndim
, nrow
, ncol
, nharm
, nvsite
;
330 int version
, generation
;
331 real value
, omega
, nu
;
332 real factor_gmx_to_omega2
;
333 real factor_omega_to_wavenumber
;
334 real
*spectrum
= NULL
;
336 gmx_output_env_t
*oenv
;
337 const char *qcleg
[] = {
338 "Heat Capacity cV (J/mol K)",
339 "Enthalpy H (kJ/mol)"
341 real
* full_hessian
= NULL
;
342 gmx_sparsematrix_t
* sparse_hessian
= NULL
;
345 { efMTX
, "-f", "hessian", ffREAD
},
346 { efTPR
, NULL
, NULL
, ffREAD
},
347 { efXVG
, "-of", "eigenfreq", ffWRITE
},
348 { efXVG
, "-ol", "eigenval", ffWRITE
},
349 { efXVG
, "-os", "spectrum", ffOPTWR
},
350 { efXVG
, "-qc", "quant_corr", ffOPTWR
},
351 { efTRN
, "-v", "eigenvec", ffWRITE
}
353 #define NFILE asize(fnm)
355 if (!parse_common_args(&argc
, argv
, 0,
356 NFILE
, fnm
, asize(pa
), pa
, asize(desc
), desc
, 0, NULL
, &oenv
))
361 /* Read tpr file for volume and number of harmonic terms */
362 read_tpxheader(ftp2fn(efTPR
, NFILE
, fnm
), &tpx
, TRUE
, &version
, &generation
);
363 snew(top_x
, tpx
.natoms
);
365 read_tpx(ftp2fn(efTPR
, NFILE
, fnm
), NULL
, box
, &natoms
,
369 nharm
= get_nharm(&mtop
, &nvsite
);
376 top
= gmx_mtop_t_to_t_topology(&mtop
);
381 if (opt2bSet("-qc", NFILE
, fnm
))
383 begin
= 7+DIM
*nvsite
;
394 printf("Using begin = %d and end = %d\n", begin
, end
);
396 /*open Hessian matrix */
397 gmx_mtxio_read(ftp2fn(efMTX
, NFILE
, fnm
), &nrow
, &ncol
, &full_hessian
, &sparse_hessian
);
399 /* Memory for eigenvalues and eigenvectors (begin..end) */
400 snew(eigenvalues
, nrow
);
401 snew(eigenvectors
, nrow
*(end
-begin
+1));
403 /* If the Hessian is in sparse format we can calculate max (ndim-1) eigenvectors,
404 * and they must start at the first one. If this is not valid we convert to full matrix
405 * storage, but warn the user that we might run out of memory...
407 if ((sparse_hessian
!= NULL
) && (begin
!= 1 || end
== ndim
))
411 fprintf(stderr
, "Cannot use sparse Hessian with first eigenvector != 1.\n");
413 else if (end
== ndim
)
415 fprintf(stderr
, "Cannot use sparse Hessian to calculate all eigenvectors.\n");
418 fprintf(stderr
, "Will try to allocate memory and convert to full matrix representation...\n");
420 snew(full_hessian
, nrow
*ncol
);
421 for (i
= 0; i
< nrow
*ncol
; i
++)
426 for (i
= 0; i
< sparse_hessian
->nrow
; i
++)
428 for (j
= 0; j
< sparse_hessian
->ndata
[i
]; j
++)
430 k
= sparse_hessian
->data
[i
][j
].col
;
431 value
= sparse_hessian
->data
[i
][j
].value
;
432 full_hessian
[i
*ndim
+k
] = value
;
433 full_hessian
[k
*ndim
+i
] = value
;
436 gmx_sparsematrix_destroy(sparse_hessian
);
437 sparse_hessian
= NULL
;
438 fprintf(stderr
, "Converted sparse to full matrix storage.\n");
441 if (full_hessian
!= NULL
)
443 /* Using full matrix storage */
444 nma_full_hessian(full_hessian
, nrow
, bM
, &top
, begin
, end
,
445 eigenvalues
, eigenvectors
);
449 /* Sparse memory storage, allocate memory for eigenvectors */
450 snew(eigenvectors
, ncol
*end
);
451 nma_sparse_hessian(sparse_hessian
, bM
, &top
, end
, eigenvalues
, eigenvectors
);
454 /* check the output, first 6 eigenvalues should be reasonably small */
456 for (i
= begin
-1; (i
< 6); i
++)
458 if (std::abs(eigenvalues
[i
]) > 1.0e-3)
465 fprintf(stderr
, "\nOne of the lowest 6 eigenvalues has a non-zero value.\n");
466 fprintf(stderr
, "This could mean that the reference structure was not\n");
467 fprintf(stderr
, "properly energy minimized.\n");
470 /* now write the output */
471 fprintf (stderr
, "Writing eigenvalues...\n");
472 out
= xvgropen(opt2fn("-ol", NFILE
, fnm
),
473 "Eigenvalues", "Eigenvalue index", "Eigenvalue [Gromacs units]",
475 if (output_env_get_print_xvgr_codes(oenv
))
479 fprintf(out
, "@ subtitle \"mass weighted\"\n");
483 fprintf(out
, "@ subtitle \"not mass weighted\"\n");
487 for (i
= 0; i
<= (end
-begin
); i
++)
489 fprintf (out
, "%6d %15g\n", begin
+i
, eigenvalues
[i
]);
494 if (opt2bSet("-qc", NFILE
, fnm
))
496 qc
= xvgropen(opt2fn("-qc", NFILE
, fnm
), "Quantum Corrections", "Eigenvector index", "", oenv
);
497 xvgr_legend(qc
, asize(qcleg
), qcleg
, oenv
);
504 printf("Writing eigenfrequencies - negative eigenvalues will be set to zero.\n");
506 out
= xvgropen(opt2fn("-of", NFILE
, fnm
),
507 "Eigenfrequencies", "Eigenvector index", "Wavenumber [cm\\S-1\\N]",
509 if (output_env_get_print_xvgr_codes(oenv
))
513 fprintf(out
, "@ subtitle \"mass weighted\"\n");
517 fprintf(out
, "@ subtitle \"not mass weighted\"\n");
522 if (opt2bSet("-os", NFILE
, fnm
) && (maxspec
> 0))
524 snew(spectrum
, maxspec
);
525 spec
= xvgropen(opt2fn("-os", NFILE
, fnm
),
526 "Vibrational spectrum based on harmonic approximation",
527 "\\f{12}w\\f{4} (cm\\S-1\\N)",
528 "Intensity [Gromacs units]",
530 for (i
= 0; (i
< maxspec
); i
++)
536 /* Gromacs units are kJ/(mol*nm*nm*amu),
537 * where amu is the atomic mass unit.
539 * For the eigenfrequencies we want to convert this to spectroscopic absorption
540 * wavenumbers given in cm^(-1), which is the frequency divided by the speed of
541 * light. Do this by first converting to omega^2 (units 1/s), take the square
542 * root, and finally divide by the speed of light (nm/ps in gromacs).
544 factor_gmx_to_omega2
= 1.0E21
/(AVOGADRO
*AMU
);
545 factor_omega_to_wavenumber
= 1.0E-5/(2.0*M_PI
*SPEED_OF_LIGHT
);
547 for (i
= begin
; (i
<= end
); i
++)
549 value
= eigenvalues
[i
-begin
];
554 omega
= std::sqrt(value
*factor_gmx_to_omega2
);
555 nu
= 1e-12*omega
/(2*M_PI
);
556 value
= omega
*factor_omega_to_wavenumber
;
557 fprintf (out
, "%6d %15g\n", i
, value
);
560 wfac
= eigenvalues
[i
-begin
]/(width
*std::sqrt(2*M_PI
));
561 for (j
= 0; (j
< maxspec
); j
++)
563 spectrum
[j
] += wfac
*std::exp(-sqr(j
-value
)/(2*sqr(width
)));
568 qcv
= cv_corr(nu
, T
);
575 fprintf (qc
, "%6d %15g %15g\n", i
, qcv
, qu
);
583 for (j
= 0; (j
< maxspec
); j
++)
585 fprintf(spec
, "%10g %10g\n", 1.0*j
, spectrum
[j
]);
591 printf("Quantum corrections for harmonic degrees of freedom\n");
592 printf("Use appropriate -first and -last options to get reliable results.\n");
593 printf("There were %d constraints and %d vsites in the simulation\n",
595 printf("Total correction to cV = %g J/mol K\n", qcvtot
);
596 printf("Total correction to H = %g kJ/mol\n", qutot
);
598 please_cite(stdout
, "Caleman2011b");
600 /* Writing eigenvectors. Note that if mass scaling was used, the eigenvectors
601 * were scaled back from mass weighted cartesian to plain cartesian in the
602 * nma_full_hessian() or nma_sparse_hessian() routines. Mass scaled vectors
603 * will not be strictly orthogonal in plain cartesian scalar products.
605 write_eigenvectors(opt2fn("-v", NFILE
, fnm
), natoms
, eigenvectors
, FALSE
, begin
, end
,
606 eWXR_NO
, NULL
, FALSE
, top_x
, bM
, eigenvalues
);