9fcaed86dc8ac4bf0d2f8a214d4caa1e547c08d1
[alexxy/gromacs.git] / src / gromacs / mdlib / qm_gamess.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  * 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #ifdef GMX_QMMM_GAMESS
42
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
47
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/legacyheaders/force.h"
54 #include "gromacs/fileio/confio.h"
55 #include "gromacs/legacyheaders/names.h"
56 #include "gromacs/legacyheaders/network.h"
57 #include "gromacs/legacyheaders/ns.h"
58 #include "gromacs/legacyheaders/nrnb.h"
59 #include "gromacs/legacyheaders/txtdump.h"
60 #include "gromacs/legacyheaders/qmmm.h"
61 #include "gromacs/utility/fatalerror.h"
62
63
64 /* QMMM sub routines */
65 /* mopac interface routines */
66
67
68 void
69 F77_FUNC(inigms, IMIGMS) (void);
70
71 void
72 F77_FUNC(endgms, ENDGMS) (void);
73
74 void
75 F77_FUNC(grads, GRADS) (int *nrqmat, real *qmcrd, int *nrmmat, real *mmchrg,
76                         real *mmcrd, real *qmgrad, real *mmgrad, real *energy);
77
78
79
80 void init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
81 {
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)
89      */
90     int
91         i, j, rank;
92     FILE
93        *out;
94     char
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"
101     };
102
103     if (PAR(cr))
104     {
105
106         if (MASTER(cr))
107         {
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.
111              */
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++)
116             {
117 #ifdef DOUBLE
118                 fprintf(out, "%10.7lf  %10.7lf  %10.7lf  %5.3lf  %2s\n",
119                         i/2.,
120                         i/3.,
121                         i/4.,
122                         qm->atomicnumberQM[i]*1.0,
123                         periodic_system[qm->atomicnumberQM[i]]);
124 #else
125                 fprintf(out, "%10.7f  %10.7f  %10.7f  %5.3f  %2s\n",
126                         i/2.,
127                         i/3.,
128                         i/4.,
129                         qm->atomicnumberQM[i]*1.0,
130                         periodic_system[qm->atomicnumberQM[i]]);
131 #endif
132             }
133             if (mm->nrMMatoms)
134             {
135                 for (j = i; j < i+2; j++)
136                 {
137 #ifdef DOUBLE
138                     fprintf(out, "%10.7lf  %10.7lf  %10.7lf  %5.3lf  BQ\n",
139                             j/5.,
140                             j/6.,
141                             j/7.,
142                             1.0);
143 #else
144                     fprintf(out, "%10.7f  %10.7f  %10.7f  %5.3f  BQ\n",
145                             j/5.,
146                             j/6.,
147                             j/7.,
148                             2.0);
149 #endif
150                 }
151             }
152             if (!qm->bTS)
153             {
154                 fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
155                         eQMbasis_names[qm->QMbasis],
156                         eQMmethod_names[qm->QMmethod]); /* see enum.h */
157             }
158             else
159             {
160                 fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
161                         eQMbasis_names[qm->QMbasis],
162                         eQMmethod_names[qm->QMmethod]); /* see enum.h */
163             }
164             fclose(out);
165         }
166         gmx_barrier(cr);
167         F77_FUNC(inigms, IMIGMS) ();
168     }
169     else /* normal serial run */
170
171     {
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.
175          */
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++)
179         {
180 #ifdef DOUBLE
181             fprintf(out, "%10.7lf  %10.7lf  %10.7lf  %5.3lf  %2s\n",
182                     i/2.,
183                     i/3.,
184                     i/4.,
185                     qm->atomicnumberQM[i]*1.0,
186                     periodic_system[qm->atomicnumberQM[i]]);
187 #else
188             fprintf(out, "%10.7f  %10.7f  %10.7f  %5.3f  %2s\n",
189                     i/2.,
190                     i/3.,
191                     i/4.,
192                     qm->atomicnumberQM[i]*1.0,
193                     periodic_system[qm->atomicnumberQM[i]]);
194 #endif
195         }
196         if (mm->nrMMatoms)
197         {
198             for (j = i; j < i+2; j++)
199             {
200 #ifdef DOUBLE
201                 fprintf(out, "%10.7lf  %10.7lf  %10.7lf  %5.3lf  BQ\n",
202                         j/5.,
203                         j/6.,
204                         j/7.,
205                         1.0);
206 #else
207                 fprintf(out, "%10.7f  %10.7f  %10.7f  %5.3f  BQ\n",
208                         j/5.,
209                         j/6.,
210                         j/7.,
211                         2.0);
212 #endif
213             }
214         }
215         if (!qm->bTS)
216         {
217             fprintf(out, "END\nBASIS %s\nRUNTYPE GRADIENT\nSCFTYPE %s\n",
218                     eQMbasis_names[qm->QMbasis],
219                     eQMmethod_names[qm->QMmethod]); /* see enum.h */
220         }
221         else
222         {
223             fprintf(out, "END\nBASIS %s\nRUNTYPE SADDLE\nSCFTYPE %s\n",
224                     eQMbasis_names[qm->QMbasis],
225                     eQMmethod_names[qm->QMmethod]); /* see enum.h */
226         }
227         F77_FUNC(inigms, IMIGMS) ();
228     }
229 }
230
231 real call_gamess(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
232                  rvec f[], rvec fshift[])
233 {
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
238      */
239     int
240         i, j, rank;
241     real
242         QMener = 0.0, *qmgrad, *mmgrad, *mmcrd, *qmcrd, energy;
243     t_QMMMrec
244        *qr;
245
246     /* copy the QMMMrec pointer */
247     qr = fr->qr;
248     snew(qmcrd, 3*(qm->nrQMatoms));
249     snew(mmcrd, 3*(mm->nrMMatoms));
250     snew(qmgrad, 3*(qm->nrQMatoms));
251     snew(mmgrad, 3*(mm->nrMMatoms));
252
253     /* copy the data from qr into the arrays that are going to be used
254      * in the fortran routines of gamess
255      */
256     for (i = 0; i < qm->nrQMatoms; i++)
257     {
258         for (j = 0; j < DIM; j++)
259         {
260             qmcrd[DIM*i+j] = 1/BOHR2NM*qm->xQM[i][j];
261         }
262     }
263     for (i = 0; i < mm->nrMMatoms; i++)
264     {
265         for (j = 0; j < DIM; j++)
266         {
267             mmcrd[DIM*i+j] = 1/BOHR2NM*mm->xMM[i][j];
268         }
269     }
270     for (i = 0; i < 3*qm->nrQMatoms; i += 3)
271     {
272         fprintf(stderr, "%8.5f, %8.5f, %8.5f\n",
273                 qmcrd[i],
274                 qmcrd[i+1],
275                 qmcrd[i+2]);
276     }
277
278     F77_FUNC(grads, GRADS) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mm->MMcharges,
279                             mmcrd, qmgrad, mmgrad, &energy);
280
281     for (i = 0; i < qm->nrQMatoms; i++)
282     {
283         for (j = 0; j < DIM; j++)
284         {
285             f[i][j]      = HARTREE_BOHR2MD*qmgrad[3*i+j];
286             fshift[i][j] = HARTREE_BOHR2MD*qmgrad[3*i+j];
287         }
288     }
289     for (i = 0; i < mm->nrMMatoms; i++)
290     {
291         for (j = 0; j < DIM; j++)
292         {
293             f[i][j]      = HARTREE_BOHR2MD*mmgrad[3*i+j];
294             fshift[i][j] = HARTREE_BOHR2MD*mmgrad[3*i+j];
295         }
296     }
297     /* convert a.u to kJ/mol */
298     QMener = energy*HARTREE2KJ*AVOGADRO;
299     return(QMener);
300 }
301
302 #else
303 int
304     gmx_qmmm_gamess_empty;
305 #endif