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
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
55 #include "gromacs/fileio/confio.h"
63 #include "gromacs/utility/fatalerror.h"
66 /* QMMM sub routines */
67 /* mopac interface routines */
71 F77_FUNC(inigms, IMIGMS) (void);
74 F77_FUNC(endgms, ENDGMS) (void);
77 F77_FUNC(grads, GRADS) (int *nrqmat, real *qmcrd, int *nrmmat, real *mmchrg,
78 real *mmcrd, real *qmgrad, real *mmgrad, real *energy);
82 void init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
84 /* it works hopelessly complicated :-)
85 * first a file is written. Then the standard gamess input/output
86 * routine is called (no system()!) to set up all fortran arrays.
87 * this routine writes a punch file, like in a normal gamess run.
88 * via this punch file the other games routines, needed for gradient
89 * and energy evaluations are called. This setup works fine for
90 * dynamics simulations. 7-6-2002 (London)
97 periodic_system[37][3] = {
98 "XX", "H ", "He", "Li", "Be", "B ", "C ", "N ",
99 "O ", "F ", "Ne", "Na", "Mg", "Al", "Si", "P ",
100 "S ", "Cl", "Ar", "K ", "Ca", "Sc", "Ti", "V ",
101 "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn", "Ga",
102 "Ge", "As", "Se", "Br", "Kr"
110 out = fopen("FOR009", "w");
111 /* of these options I am not completely sure.... the overall
112 * preformance on more than 4 cpu's is rather poor at the moment.
114 fprintf(out, "memory 48000000\nPARALLEL IOMODE SCREENED\n");
115 fprintf(out, "ELEC %d\nMULT %d\nSUPER ON\nNOSYM\nGEOMETRY ANGSTROM\n",
116 qm->nelectrons, qm->multiplicity);
117 for (i = 0; i < qm->nrQMatoms; i++)
120 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf %2s\n",
124 qm->atomicnumberQM[i]*1.0,
125 periodic_system[qm->atomicnumberQM[i]]);
127 fprintf(out, "%10.7f %10.7f %10.7f %5.3f %2s\n",
131 qm->atomicnumberQM[i]*1.0,
132 periodic_system[qm->atomicnumberQM[i]]);
137 for (j = i; j < i+2; j++)
140 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf BQ\n",
146 fprintf(out, "%10.7f %10.7f %10.7f %5.3f BQ\n",
156 fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
157 eQMbasis_names[qm->QMbasis],
158 eQMmethod_names[qm->QMmethod]); /* see enum.h */
162 fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
163 eQMbasis_names[qm->QMbasis],
164 eQMmethod_names[qm->QMmethod]); /* see enum.h */
169 F77_FUNC(inigms, IMIGMS) ();
171 else /* normal serial run */
174 out = fopen("FOR009", "w");
175 /* of these options I am not completely sure.... the overall
176 * preformance on more than 4 cpu's is rather poor at the moment.
178 fprintf(out, "ELEC %d\nMULT %d\nSUPER ON\nNOSYM\nGEOMETRY ANGSTROM\n",
179 qm->nelectrons, qm->multiplicity);
180 for (i = 0; i < qm->nrQMatoms; i++)
183 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf %2s\n",
187 qm->atomicnumberQM[i]*1.0,
188 periodic_system[qm->atomicnumberQM[i]]);
190 fprintf(out, "%10.7f %10.7f %10.7f %5.3f %2s\n",
194 qm->atomicnumberQM[i]*1.0,
195 periodic_system[qm->atomicnumberQM[i]]);
200 for (j = i; j < i+2; j++)
203 fprintf(out, "%10.7lf %10.7lf %10.7lf %5.3lf BQ\n",
209 fprintf(out, "%10.7f %10.7f %10.7f %5.3f BQ\n",
219 fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
220 eQMbasis_names[qm->QMbasis],
221 eQMmethod_names[qm->QMmethod]); /* see enum.h */
225 fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
226 eQMbasis_names[qm->QMbasis],
227 eQMmethod_names[qm->QMmethod]); /* see enum.h */
229 F77_FUNC(inigms, IMIGMS) ();
233 real call_gamess(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
234 rvec f[], rvec fshift[])
236 /* do the actual QMMM calculation using GAMESS-UK. In this
237 * implementation (3-2001) a system call is made to the GAMESS-UK
238 * binary. Now we are working to get the electron integral, SCF, and
239 * gradient routines linked directly
244 QMener = 0.0, *qmgrad, *mmgrad, *mmcrd, *qmcrd, energy;
248 /* copy the QMMMrec pointer */
250 snew(qmcrd, 3*(qm->nrQMatoms));
251 snew(mmcrd, 3*(mm->nrMMatoms));
252 snew(qmgrad, 3*(qm->nrQMatoms));
253 snew(mmgrad, 3*(mm->nrMMatoms));
255 /* copy the data from qr into the arrays that are going to be used
256 * in the fortran routines of gamess
258 for (i = 0; i < qm->nrQMatoms; i++)
260 for (j = 0; j < DIM; j++)
262 qmcrd[DIM*i+j] = 1/BOHR2NM*qm->xQM[i][j];
265 for (i = 0; i < mm->nrMMatoms; i++)
267 for (j = 0; j < DIM; j++)
269 mmcrd[DIM*i+j] = 1/BOHR2NM*mm->xMM[i][j];
272 for (i = 0; i < 3*qm->nrQMatoms; i += 3)
274 fprintf(stderr, "%8.5f, %8.5f, %8.5f\n",
280 F77_FUNC(grads, GRADS) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mm->MMcharges,
281 mmcrd, qmgrad, mmgrad, &energy);
283 for (i = 0; i < qm->nrQMatoms; i++)
285 for (j = 0; j < DIM; j++)
287 f[i][j] = HARTREE_BOHR2MD*qmgrad[3*i+j];
288 fshift[i][j] = HARTREE_BOHR2MD*qmgrad[3*i+j];
291 for (i = 0; i < mm->nrMMatoms; i++)
293 for (j = 0; j < DIM; j++)
295 f[i][j] = HARTREE_BOHR2MD*mmgrad[3*i+j];
296 fshift[i][j] = HARTREE_BOHR2MD*mmgrad[3*i+j];
299 /* convert a.u to kJ/mol */
300 QMener = energy*HARTREE2KJ*AVOGADRO;
306 gmx_qmmm_gamess_empty;