Merge branch 'release-4-6', adds the nbnxn functionality
[alexxy/gromacs.git] / src / gromacs / mdlib / qm_mopac.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14  
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROwing Monsters And Cloning Shrimps
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #ifdef GMX_QMMM_MOPAC
40
41 #include <math.h>
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "macros.h"
45 #include "smalloc.h"
46 #include "physics.h"
47 #include "macros.h"
48 #include "vec.h"
49 #include "force.h"
50 #include "invblock.h"
51 #include "confio.h"
52 #include "names.h"
53 #include "network.h"
54 #include "pbc.h"
55 #include "ns.h"
56 #include "nrnb.h"
57 #include "bondf.h"
58 #include "mshift.h"
59 #include "txtdump.h"
60 #include "copyrite.h"
61 #include "qmmm.h"
62 #include <stdio.h>
63 #include <string.h>
64 #include "gmx_fatal.h"
65 #include "typedefs.h"
66 #include <stdlib.h>
67
68
69 /* mopac interface routines */
70
71 #ifndef F77_FUNC
72 #define F77_FUNC(name,NAME) name ## _
73 #endif
74
75
76 void 
77 F77_FUNC(domldt,DOMLDT)(int *nrqmat, int labels[], char keywords[]);
78
79 void 
80 F77_FUNC(domop,DOMOP)(int *nrqmat,double qmcrd[],int *nrmmat,
81                       double mmchrg[],double mmcrd[],double qmgrad[],
82                       double mmgrad[], double *energy,double qmcharges[]);
83
84
85
86 void init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
87 {
88   /* initializes the mopac routines ans sets up the semiempirical
89    * computation by calling moldat(). The inline mopac routines can
90    * only perform gradient operations. If one would like to optimize a
91    * structure or find a transition state at PM3 level, gaussian is
92    * used instead.
93    */
94   char 
95     *keywords;
96   
97   snew(keywords,240);
98   
99   if(!qm->bSH){    /* if rerun then grad should not be done! */
100     sprintf(keywords,"PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
101             qm->QMcharge,
102             eQMmethod_names[qm->QMmethod]);
103   }
104   else
105     sprintf(keywords,"PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
106             qm->QMcharge,
107             eQMmethod_names[qm->QMmethod],
108             qm->CASorbitals,qm->CASelectrons/2);
109   F77_FUNC(domldt,DOMLDT)(&qm->nrQMatoms,qm->atomicnumberQM,keywords);
110   fprintf(stderr,"keywords are: %s\n",keywords);
111   free(keywords);
112   
113 } /* init_mopac */
114
115 real call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, 
116                 rvec f[], rvec fshift[])
117 {
118   /* do the actual QMMM calculation using directly linked mopac subroutines 
119    */
120   double /* always double as the MOPAC routines are always compiled in
121             double precission! */
122     *qmcrd=NULL,*qmchrg=NULL,*mmcrd=NULL,*mmchrg=NULL,
123     *qmgrad,*mmgrad=NULL,energy;
124   int
125     i,j;
126   real
127     QMener=0.0;
128   snew(qmcrd, 3*(qm->nrQMatoms));
129   snew(qmgrad,3*(qm->nrQMatoms));
130   /* copy the data from qr into the arrays that are going to be used
131    * in the fortran routines of MOPAC
132    */
133   for(i=0;i<qm->nrQMatoms;i++){
134     for (j=0;j<DIM;j++){
135       qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
136     }
137   }
138   if(mm->nrMMatoms){
139     /* later we will add the point charges here. There are some
140      * conceptual problems with semi-empirical QM in combination with
141      * point charges that we need to solve first....  
142      */
143     gmx_fatal(FARGS,"At present only ONIOM is allowed in combination"
144                 " with MOPAC QM subroutines\n");
145   }
146   else {
147     /* now compute the energy and the gradients.
148      */
149       
150     snew(qmchrg,qm->nrQMatoms);    
151     F77_FUNC(domop,DOMOP)(&qm->nrQMatoms,qmcrd,&mm->nrMMatoms,
152            mmchrg,mmcrd,qmgrad,mmgrad,&energy,qmchrg);
153     /* add the gradients to the f[] array, and also to the fshift[].
154      * the mopac gradients are in kCal/angstrom.
155      */
156     for(i=0;i<qm->nrQMatoms;i++){
157       for(j=0;j<DIM;j++){
158         f[i][j]       = (real)10*CAL2JOULE*qmgrad[3*i+j];
159         fshift[i][j]  = (real)10*CAL2JOULE*qmgrad[3*i+j];
160       }
161     }
162     QMener = (real)CAL2JOULE*energy;
163     /* do we do something with the mulliken charges?? */
164
165     free(qmchrg);
166 }
167   free(qmgrad);
168   free(qmcrd);
169   return (QMener);
170 }
171
172 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, 
173                    rvec f[], rvec fshift[])
174 {
175   /* do the actual SH QMMM calculation using directly linked mopac
176    subroutines */
177
178   double /* always double as the MOPAC routines are always compiled in
179             double precission! */
180     *qmcrd=NULL,*qmchrg=NULL,*mmcrd=NULL,*mmchrg=NULL,
181     *qmgrad,*mmgrad=NULL,energy;
182   int
183     i,j;
184   real
185     QMener=0.0;
186
187   snew(qmcrd, 3*(qm->nrQMatoms));
188   snew(qmgrad,3*(qm->nrQMatoms));
189   /* copy the data from qr into the arrays that are going to be used
190    * in the fortran routines of MOPAC
191    */
192   for(i=0;i<qm->nrQMatoms;i++){
193     for (j=0;j<DIM;j++){
194       qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
195     }
196   }
197   if(mm->nrMMatoms){
198     /* later we will add the point charges here. There are some
199      * conceptual problems with semi-empirical QM in combination with
200      * point charges that we need to solve first....  
201      */
202     gmx_fatal(FARGS,"At present only ONIOM is allowed in combination with MOPAC\n");
203   }
204   else {
205     /* now compute the energy and the gradients.
206      */
207     snew(qmchrg,qm->nrQMatoms);    
208
209     F77_FUNC(domop,DOMOP)(&qm->nrQMatoms,qmcrd,&mm->nrMMatoms,
210            mmchrg,mmcrd,qmgrad,mmgrad,&energy,qmchrg);
211     /* add the gradients to the f[] array, and also to the fshift[].
212      * the mopac gradients are in kCal/angstrom.
213      */
214     for(i=0;i<qm->nrQMatoms;i++){
215       for(j=0;j<DIM;j++){
216         f[i][j]      = (real)10*CAL2JOULE*qmgrad[3*i+j];
217         fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
218       }
219     }
220     QMener = (real)CAL2JOULE*energy;
221   }
222   free(qmgrad);
223   free(qmcrd);
224   return (QMener);
225 } /* call_mopac_SH */
226
227 #else
228 int
229 gmx_qmmm_mopac_empty;
230 #endif