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, 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.
41 #ifdef GMX_QMMM_GAMESS
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/nrnb.h"
54 #include "gromacs/legacyheaders/ns.h"
55 #include "gromacs/legacyheaders/qmmm.h"
56 #include "gromacs/legacyheaders/txtdump.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
64 /* QMMM sub routines */
65 /* mopac interface routines */
69 F77_FUNC(inigms, IMIGMS) (void);
72 F77_FUNC(endgms, ENDGMS) (void);
75 F77_FUNC(grads, GRADS) (int *nrqmat, real *qmcrd, int *nrmmat, real *mmchrg,
76 real *mmcrd, real *qmgrad, real *mmgrad, real *energy);
80 void init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
82 /* it works hopelessly complicated :-)
83 * first a file is written. Then the standard gamess input/output
84 * routine is called (no system()!) to set up all fortran arrays.
85 * this routine writes a punch file, like in a normal gamess run.
86 * via this punch file the other games routines, needed for gradient
87 * and energy evaluations are called. This setup works fine for
88 * dynamics simulations. 7-6-2002 (London)
95 periodic_system[37][3] = {
96 "XX", "H ", "He", "Li", "Be", "B ", "C ", "N ",
97 "O ", "F ", "Ne", "Na", "Mg", "Al", "Si", "P ",
98 "S ", "Cl", "Ar", "K ", "Ca", "Sc", "Ti", "V ",
99 "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn", "Ga",
100 "Ge", "As", "Se", "Br", "Kr"
108 out = fopen("FOR009", "w");
109 /* of these options I am not completely sure.... the overall
110 * preformance on more than 4 cpu's is rather poor at the moment.
112 fprintf(out, "memory 48000000\nPARALLEL IOMODE SCREENED\n");
113 fprintf(out, "ELEC %d\nMULT %d\nSUPER ON\nNOSYM\nGEOMETRY ANGSTROM\n",
114 qm->nelectrons, qm->multiplicity);
115 for (i = 0; i < qm->nrQMatoms; i++)
118 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf %2s\n",
122 qm->atomicnumberQM[i]*1.0,
123 periodic_system[qm->atomicnumberQM[i]]);
125 fprintf(out, "%10.7f %10.7f %10.7f %5.3f %2s\n",
129 qm->atomicnumberQM[i]*1.0,
130 periodic_system[qm->atomicnumberQM[i]]);
135 for (j = i; j < i+2; j++)
138 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf BQ\n",
144 fprintf(out, "%10.7f %10.7f %10.7f %5.3f BQ\n",
154 fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
155 eQMbasis_names[qm->QMbasis],
156 eQMmethod_names[qm->QMmethod]); /* see enum.h */
160 fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
161 eQMbasis_names[qm->QMbasis],
162 eQMmethod_names[qm->QMmethod]); /* see enum.h */
167 F77_FUNC(inigms, IMIGMS) ();
169 else /* normal serial run */
172 out = fopen("FOR009", "w");
173 /* of these options I am not completely sure.... the overall
174 * preformance on more than 4 cpu's is rather poor at the moment.
176 fprintf(out, "ELEC %d\nMULT %d\nSUPER ON\nNOSYM\nGEOMETRY ANGSTROM\n",
177 qm->nelectrons, qm->multiplicity);
178 for (i = 0; i < qm->nrQMatoms; i++)
181 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf %2s\n",
185 qm->atomicnumberQM[i]*1.0,
186 periodic_system[qm->atomicnumberQM[i]]);
188 fprintf(out, "%10.7f %10.7f %10.7f %5.3f %2s\n",
192 qm->atomicnumberQM[i]*1.0,
193 periodic_system[qm->atomicnumberQM[i]]);
198 for (j = i; j < i+2; j++)
201 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf BQ\n",
207 fprintf(out, "%10.7f %10.7f %10.7f %5.3f BQ\n",
217 fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
218 eQMbasis_names[qm->QMbasis],
219 eQMmethod_names[qm->QMmethod]); /* see enum.h */
223 fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
224 eQMbasis_names[qm->QMbasis],
225 eQMmethod_names[qm->QMmethod]); /* see enum.h */
227 F77_FUNC(inigms, IMIGMS) ();
231 real call_gamess(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
232 rvec f[], rvec fshift[])
234 /* do the actual QMMM calculation using GAMESS-UK. In this
235 * implementation (3-2001) a system call is made to the GAMESS-UK
236 * binary. Now we are working to get the electron integral, SCF, and
237 * gradient routines linked directly
242 QMener = 0.0, *qmgrad, *mmgrad, *mmcrd, *qmcrd, energy;
246 /* copy the QMMMrec pointer */
248 snew(qmcrd, 3*(qm->nrQMatoms));
249 snew(mmcrd, 3*(mm->nrMMatoms));
250 snew(qmgrad, 3*(qm->nrQMatoms));
251 snew(mmgrad, 3*(mm->nrMMatoms));
253 /* copy the data from qr into the arrays that are going to be used
254 * in the fortran routines of gamess
256 for (i = 0; i < qm->nrQMatoms; i++)
258 for (j = 0; j < DIM; j++)
260 qmcrd[DIM*i+j] = 1/BOHR2NM*qm->xQM[i][j];
263 for (i = 0; i < mm->nrMMatoms; i++)
265 for (j = 0; j < DIM; j++)
267 mmcrd[DIM*i+j] = 1/BOHR2NM*mm->xMM[i][j];
270 for (i = 0; i < 3*qm->nrQMatoms; i += 3)
272 fprintf(stderr, "%8.5f, %8.5f, %8.5f\n",
278 F77_FUNC(grads, GRADS) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mm->MMcharges,
279 mmcrd, qmgrad, mmgrad, &energy);
281 for (i = 0; i < qm->nrQMatoms; i++)
283 for (j = 0; j < DIM; j++)
285 f[i][j] = HARTREE_BOHR2MD*qmgrad[3*i+j];
286 fshift[i][j] = HARTREE_BOHR2MD*qmgrad[3*i+j];
289 for (i = 0; i < mm->nrMMatoms; i++)
291 for (j = 0; j < DIM; j++)
293 f[i][j] = HARTREE_BOHR2MD*mmgrad[3*i+j];
294 fshift[i][j] = HARTREE_BOHR2MD*mmgrad[3*i+j];
297 /* convert a.u to kJ/mol */
298 QMener = energy*HARTREE2KJ*AVOGADRO;
304 gmx_qmmm_gamess_empty;