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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * 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.
68 #include "gmx_fatal.h"
73 /* mopac interface routines */
75 F77_FUNC(domldt,DOMLDT)(int *nrqmat, int labels[], char keywords[]);
78 F77_FUNC(domop,DOMOP)(int *nrqmat,double qmcrd[],int *nrmmat,
79 double mmchrg[],double mmcrd[],double qmgrad[],
80 double mmgrad[], double *energy,double qmcharges[]);
84 void init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
86 /* initializes the mopac routines ans sets up the semiempirical
87 * computation by calling moldat(). The inline mopac routines can
88 * only perform gradient operations. If one would like to optimize a
89 * structure or find a transition state at PM3 level, gaussian is
97 if(!qm->bSH){ /* if rerun then grad should not be done! */
98 sprintf(keywords,"PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
100 eQMmethod_names[qm->QMmethod]);
103 sprintf(keywords,"PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
105 eQMmethod_names[qm->QMmethod],
106 qm->CASorbitals,qm->CASelectrons/2);
107 F77_FUNC(domldt,DOMLDT)(&qm->nrQMatoms,qm->atomicnumberQM,keywords);
108 fprintf(stderr,"keywords are: %s\n",keywords);
113 real call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
114 rvec f[], rvec fshift[])
116 /* do the actual QMMM calculation using directly linked mopac subroutines
118 double /* always double as the MOPAC routines are always compiled in
119 double precission! */
120 *qmcrd=NULL,*qmchrg=NULL,*mmcrd=NULL,*mmchrg=NULL,
121 *qmgrad,*mmgrad=NULL,energy;
126 snew(qmcrd, 3*(qm->nrQMatoms));
127 snew(qmgrad,3*(qm->nrQMatoms));
128 /* copy the data from qr into the arrays that are going to be used
129 * in the fortran routines of MOPAC
131 for(i=0;i<qm->nrQMatoms;i++){
133 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
137 /* later we will add the point charges here. There are some
138 * conceptual problems with semi-empirical QM in combination with
139 * point charges that we need to solve first....
141 gmx_fatal(FARGS,"At present only ONIOM is allowed in combination"
142 " with MOPAC QM subroutines\n");
145 /* now compute the energy and the gradients.
148 snew(qmchrg,qm->nrQMatoms);
149 F77_FUNC(domop,DOMOP)(&qm->nrQMatoms,qmcrd,&mm->nrMMatoms,
150 mmchrg,mmcrd,qmgrad,mmgrad,&energy,qmchrg);
151 /* add the gradients to the f[] array, and also to the fshift[].
152 * the mopac gradients are in kCal/angstrom.
154 for(i=0;i<qm->nrQMatoms;i++){
156 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
157 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
160 QMener = (real)CAL2JOULE*energy;
161 /* do we do something with the mulliken charges?? */
170 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
171 rvec f[], rvec fshift[])
173 /* do the actual SH QMMM calculation using directly linked mopac
176 double /* always double as the MOPAC routines are always compiled in
177 double precission! */
178 *qmcrd=NULL,*qmchrg=NULL,*mmcrd=NULL,*mmchrg=NULL,
179 *qmgrad,*mmgrad=NULL,energy;
185 snew(qmcrd, 3*(qm->nrQMatoms));
186 snew(qmgrad,3*(qm->nrQMatoms));
187 /* copy the data from qr into the arrays that are going to be used
188 * in the fortran routines of MOPAC
190 for(i=0;i<qm->nrQMatoms;i++){
192 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
196 /* later we will add the point charges here. There are some
197 * conceptual problems with semi-empirical QM in combination with
198 * point charges that we need to solve first....
200 gmx_fatal(FARGS,"At present only ONIOM is allowed in combination with MOPAC\n");
203 /* now compute the energy and the gradients.
205 snew(qmchrg,qm->nrQMatoms);
207 F77_FUNC(domop,DOMOP)(&qm->nrQMatoms,qmcrd,&mm->nrMMatoms,
208 mmchrg,mmcrd,qmgrad,mmgrad,&energy,qmchrg);
209 /* add the gradients to the f[] array, and also to the fshift[].
210 * the mopac gradients are in kCal/angstrom.
212 for(i=0;i<qm->nrQMatoms;i++){
214 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
215 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
218 QMener = (real)CAL2JOULE*energy;
223 } /* call_mopac_SH */
227 gmx_qmmm_mopac_empty;