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