Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / mdlib / qm_mopac.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #ifdef GMX_QMMM_MOPAC
43
44 #include <math.h>
45 #include "sysstuff.h"
46 #include "typedefs.h"
47 #include "macros.h"
48 #include "smalloc.h"
49 #include "assert.h"
50 #include "physics.h"
51 #include "macros.h"
52 #include "vec.h"
53 #include "force.h"
54 #include "invblock.h"
55 #include "confio.h"
56 #include "names.h"
57 #include "network.h"
58 #include "pbc.h"
59 #include "ns.h"
60 #include "nrnb.h"
61 #include "bondf.h"
62 #include "mshift.h"
63 #include "txtdump.h"
64 #include "copyrite.h"
65 #include "qmmm.h"
66 #include <stdio.h>
67 #include <string.h>
68 #include "gmx_fatal.h"
69 #include "typedefs.h"
70 #include <stdlib.h>
71
72
73 /* mopac interface routines */
74 void 
75 F77_FUNC(domldt,DOMLDT)(int *nrqmat, int labels[], char keywords[]);
76
77 void 
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[]);
81
82
83
84 void init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
85 {
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
90    * used instead.
91    */
92   char 
93     *keywords;
94   
95   snew(keywords,240);
96   
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",
99             qm->QMcharge,
100             eQMmethod_names[qm->QMmethod]);
101   }
102   else
103     sprintf(keywords,"PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
104             qm->QMcharge,
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);
109   free(keywords);
110   
111 } /* init_mopac */
112
113 real call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, 
114                 rvec f[], rvec fshift[])
115 {
116   /* do the actual QMMM calculation using directly linked mopac subroutines 
117    */
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;
122   int
123     i,j;
124   real
125     QMener=0.0;
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
130    */
131   for(i=0;i<qm->nrQMatoms;i++){
132     for (j=0;j<DIM;j++){
133       qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
134     }
135   }
136   if(mm->nrMMatoms){
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....  
140      */
141     gmx_fatal(FARGS,"At present only ONIOM is allowed in combination"
142                 " with MOPAC QM subroutines\n");
143   }
144   else {
145     /* now compute the energy and the gradients.
146      */
147       
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.
153      */
154     for(i=0;i<qm->nrQMatoms;i++){
155       for(j=0;j<DIM;j++){
156         f[i][j]       = (real)10*CAL2JOULE*qmgrad[3*i+j];
157         fshift[i][j]  = (real)10*CAL2JOULE*qmgrad[3*i+j];
158       }
159     }
160     QMener = (real)CAL2JOULE*energy;
161     /* do we do something with the mulliken charges?? */
162
163     free(qmchrg);
164 }
165   free(qmgrad);
166   free(qmcrd);
167   return (QMener);
168 }
169
170 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, 
171                    rvec f[], rvec fshift[])
172 {
173   /* do the actual SH QMMM calculation using directly linked mopac
174    subroutines */
175
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;
180   int
181     i,j;
182   real
183     QMener=0.0;
184
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
189    */
190   for(i=0;i<qm->nrQMatoms;i++){
191     for (j=0;j<DIM;j++){
192       qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
193     }
194   }
195   if(mm->nrMMatoms){
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....  
199      */
200     gmx_fatal(FARGS,"At present only ONIOM is allowed in combination with MOPAC\n");
201   }
202   else {
203     /* now compute the energy and the gradients.
204      */
205     snew(qmchrg,qm->nrQMatoms);    
206
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.
211      */
212     for(i=0;i<qm->nrQMatoms;i++){
213       for(j=0;j<DIM;j++){
214         f[i][j]      = (real)10*CAL2JOULE*qmgrad[3*i+j];
215         fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
216       }
217     }
218     QMener = (real)CAL2JOULE*energy;
219   }
220   free(qmgrad);
221   free(qmcrd);
222   return (QMener);
223 } /* call_mopac_SH */
224
225 #else
226 int
227 gmx_qmmm_mopac_empty;
228 #endif