3ae40c2a840f04633487f665dbe40ae60a0c0e71
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / convparm.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 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include <math.h>
41 #include <string.h>
42
43 #include "gromacs/math/units.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "topio.h"
49 #include "toputil.h"
50 #include "convparm.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gpp_atomtype.h"
53 #include "gromacs/math/utilities.h"
54
55 static int round_check(real r, int limit, int ftype, const char *name)
56 {
57     int i;
58
59     if (r >= 0)
60     {
61         i = (int)(r + 0.5);
62     }
63     else
64     {
65         i = (int)(r - 0.5);
66     }
67
68     if (r-i > 0.01 || r-i < -0.01)
69     {
70         gmx_fatal(FARGS, "A non-integer value (%f) was supplied for '%s' in %s",
71                   r, name, interaction_function[ftype].longname);
72     }
73
74     if (i < limit)
75     {
76         gmx_fatal(FARGS, "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
77                   name, interaction_function[ftype].longname, i, limit);
78     }
79
80     return i;
81 }
82
83 static void set_ljparams(int comb, double reppow, real v, real w,
84                          real *c6, real *c12)
85 {
86     if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
87     {
88         if (v >= 0)
89         {
90             *c6  = 4*w*pow(v, 6.0);
91             *c12 = 4*w*pow(v, reppow);
92         }
93         else
94         {
95             /* Interpret negative sigma as c6=0 and c12 with -sigma */
96             *c6  = 0;
97             *c12 = 4*w*pow(-v, reppow);
98         }
99     }
100     else
101     {
102         *c6  = v;
103         *c12 = w;
104     }
105 }
106
107 /* A return value of 0 means parameters were assigned successfully,
108  * returning -1 means this is an all-zero interaction that should not be added.
109  */
110 static int
111 assign_param(t_functype ftype, t_iparams *newparam,
112              real old[MAXFORCEPARAM], int comb, double reppow)
113 {
114     int      i, j;
115     real     tmp;
116     gmx_bool all_param_zero = TRUE;
117
118     /* Set to zero */
119     for (j = 0; (j < MAXFORCEPARAM); j++)
120     {
121         newparam->generic.buf[j] = 0.0;
122         /* If all parameters are zero we might not add some interaction types (selected below).
123          * We cannot apply this to ALL interactions, since many have valid reasons for having
124          * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
125          * we use it for angles and torsions that are typically generated automatically.
126          */
127         all_param_zero = (all_param_zero == TRUE) && fabs(old[j]) < GMX_REAL_MIN;
128     }
129
130     if (all_param_zero == TRUE)
131     {
132         if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS ||
133             ftype == F_PDIHS || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
134         {
135             return -1;
136         }
137     }
138
139     switch (ftype)
140     {
141         case F_G96ANGLES:
142             /* Post processing of input data: store cosine iso angle itself */
143             newparam->harmonic.rA  = cos(old[0]*DEG2RAD);
144             newparam->harmonic.krA = old[1];
145             newparam->harmonic.rB  = cos(old[2]*DEG2RAD);
146             newparam->harmonic.krB = old[3];
147             break;
148         case F_G96BONDS:
149             /* Post processing of input data: store square of length itself */
150             newparam->harmonic.rA  = sqr(old[0]);
151             newparam->harmonic.krA = old[1];
152             newparam->harmonic.rB  = sqr(old[2]);
153             newparam->harmonic.krB = old[3];
154             break;
155         case F_FENEBONDS:
156             newparam->fene.bm = old[0];
157             newparam->fene.kb = old[1];
158             break;
159         case F_RESTRBONDS:
160             newparam->restraint.lowA = old[0];
161             newparam->restraint.up1A = old[1];
162             newparam->restraint.up2A = old[2];
163             newparam->restraint.kA   = old[3];
164             newparam->restraint.lowB = old[4];
165             newparam->restraint.up1B = old[5];
166             newparam->restraint.up2B = old[6];
167             newparam->restraint.kB   = old[7];
168             break;
169         case F_TABBONDS:
170         case F_TABBONDSNC:
171         case F_TABANGLES:
172         case F_TABDIHS:
173             newparam->tab.table = round_check(old[0], 0, ftype, "table index");
174             newparam->tab.kA    = old[1];
175             newparam->tab.kB    = old[3];
176             break;
177         case F_CROSS_BOND_BONDS:
178             newparam->cross_bb.r1e = old[0];
179             newparam->cross_bb.r2e = old[1];
180             newparam->cross_bb.krr = old[2];
181             break;
182         case F_CROSS_BOND_ANGLES:
183             newparam->cross_ba.r1e = old[0];
184             newparam->cross_ba.r2e = old[1];
185             newparam->cross_ba.r3e = old[2];
186             newparam->cross_ba.krt = old[3];
187             break;
188         case F_UREY_BRADLEY:
189             newparam->u_b.thetaA  = old[0];
190             newparam->u_b.kthetaA = old[1];
191             newparam->u_b.r13A    = old[2];
192             newparam->u_b.kUBA    = old[3];
193             newparam->u_b.thetaB  = old[4];
194             newparam->u_b.kthetaB = old[5];
195             newparam->u_b.r13B    = old[6];
196             newparam->u_b.kUBB    = old[7];
197             break;
198         case F_QUARTIC_ANGLES:
199             newparam->qangle.theta = old[0];
200             for (i = 0; i < 5; i++)
201             {
202                 newparam->qangle.c[i] = old[i+1];
203             }
204             break;
205         case F_LINEAR_ANGLES:
206             newparam->linangle.aA    = old[0];
207             newparam->linangle.klinA = old[1];
208             newparam->linangle.aB    = old[2];
209             newparam->linangle.klinB = old[3];
210             break;
211         case F_BONDS:
212         case F_ANGLES:
213         case F_HARMONIC:
214         case F_IDIHS:
215             newparam->harmonic.rA  = old[0];
216             newparam->harmonic.krA = old[1];
217             newparam->harmonic.rB  = old[2];
218             newparam->harmonic.krB = old[3];
219             break;
220         case F_RESTRANGLES:
221             newparam->harmonic.rA  = old[0];
222             newparam->harmonic.krA = old[1];
223             break;
224         case F_MORSE:
225             newparam->morse.b0A    = old[0];
226             newparam->morse.cbA    = old[1];
227             newparam->morse.betaA  = old[2];
228             newparam->morse.b0B    = old[3];
229             newparam->morse.cbB    = old[4];
230             newparam->morse.betaB  = old[5];
231             break;
232         case F_CUBICBONDS:
233             newparam->cubic.b0    = old[0];
234             newparam->cubic.kb    = old[1];
235             newparam->cubic.kcub  = old[2];
236             break;
237         case F_CONNBONDS:
238             break;
239         case F_POLARIZATION:
240             newparam->polarize.alpha = old[0];
241             break;
242         case F_ANHARM_POL:
243             newparam->anharm_polarize.alpha = old[0];
244             newparam->anharm_polarize.drcut = old[1];
245             newparam->anharm_polarize.khyp  = old[2];
246             break;
247         case F_WATER_POL:
248             newparam->wpol.al_x   = old[0];
249             newparam->wpol.al_y   = old[1];
250             newparam->wpol.al_z   = old[2];
251             newparam->wpol.rOH    = old[3];
252             newparam->wpol.rHH    = old[4];
253             newparam->wpol.rOD    = old[5];
254             break;
255         case F_THOLE_POL:
256             newparam->thole.a      = old[0];
257             newparam->thole.alpha1 = old[1];
258             newparam->thole.alpha2 = old[2];
259             if ((old[1] > 0) && (old[2] > 0))
260             {
261                 newparam->thole.rfac = old[0]*pow(old[1]*old[2], -1.0/6.0);
262             }
263             else
264             {
265                 newparam->thole.rfac = 1;
266             }
267             break;
268         case F_BHAM:
269             newparam->bham.a = old[0];
270             newparam->bham.b = old[1];
271             newparam->bham.c = old[2];
272             break;
273         case F_LJ14:
274             set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
275             set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
276             break;
277         case F_LJC14_Q:
278             newparam->ljc14.fqq = old[0];
279             newparam->ljc14.qi  = old[1];
280             newparam->ljc14.qj  = old[2];
281             set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
282             break;
283         case F_LJC_PAIRS_NB:
284             newparam->ljcnb.qi = old[0];
285             newparam->ljcnb.qj = old[1];
286             set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
287             break;
288         case F_LJ:
289             set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
290             break;
291         case F_PDIHS:
292         case F_PIDIHS:
293         case F_ANGRES:
294         case F_ANGRESZ:
295             newparam->pdihs.phiA = old[0];
296             newparam->pdihs.cpA  = old[1];
297
298             /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
299              * so I have changed the lower limit to -99 /EL
300              */
301             newparam->pdihs.phiB = old[3];
302             newparam->pdihs.cpB  = old[4];
303             /* If both force constants are zero there is no interaction. Return -1 to signal
304              * this entry should NOT be added.
305              */
306             if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
307             {
308                 return -1;
309             }
310
311             newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
312
313             break;
314         case F_RESTRDIHS:
315             newparam->pdihs.phiA = old[0];
316             newparam->pdihs.cpA  = old[1];
317             break;
318         case F_POSRES:
319             newparam->posres.fcA[XX]   = old[0];
320             newparam->posres.fcA[YY]   = old[1];
321             newparam->posres.fcA[ZZ]   = old[2];
322             newparam->posres.fcB[XX]   = old[3];
323             newparam->posres.fcB[YY]   = old[4];
324             newparam->posres.fcB[ZZ]   = old[5];
325             newparam->posres.pos0A[XX] = old[6];
326             newparam->posres.pos0A[YY] = old[7];
327             newparam->posres.pos0A[ZZ] = old[8];
328             newparam->posres.pos0B[XX] = old[9];
329             newparam->posres.pos0B[YY] = old[10];
330             newparam->posres.pos0B[ZZ] = old[11];
331             break;
332         case F_FBPOSRES:
333             newparam->fbposres.geom     = round_check(old[0], 0, ftype, "geometry");
334             if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
335             {
336                 gmx_fatal(FARGS, "Invalid geometry for flat-bottomed position restraint.\n"
337                           "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
338                           newparam->fbposres.geom);
339             }
340             newparam->fbposres.r        = old[1];
341             newparam->fbposres.k        = old[2];
342             newparam->fbposres.pos0[XX] = old[3];
343             newparam->fbposres.pos0[YY] = old[4];
344             newparam->fbposres.pos0[ZZ] = old[5];
345             break;
346         case F_DISRES:
347             newparam->disres.label = round_check(old[0], 0, ftype, "label");
348             newparam->disres.type  = round_check(old[1], 1, ftype, "type'");
349             newparam->disres.low   = old[2];
350             newparam->disres.up1   = old[3];
351             newparam->disres.up2   = old[4];
352             newparam->disres.kfac  = old[5];
353             break;
354         case F_ORIRES:
355             newparam->orires.ex    = round_check(old[0], 1, ftype, "experiment") - 1;
356             newparam->orires.label = round_check(old[1], 1, ftype, "label");
357             newparam->orires.power = round_check(old[2], 0, ftype, "power");
358             newparam->orires.c     = old[3];
359             newparam->orires.obs   = old[4];
360             newparam->orires.kfac  = old[5];
361             break;
362         case F_DIHRES:
363             newparam->dihres.phiA   = old[0];
364             newparam->dihres.dphiA  = old[1];
365             newparam->dihres.kfacA  = old[2];
366             newparam->dihres.phiB   = old[3];
367             newparam->dihres.dphiB  = old[4];
368             newparam->dihres.kfacB  = old[5];
369             break;
370         case F_RBDIHS:
371             for (i = 0; (i < NR_RBDIHS); i++)
372             {
373                 newparam->rbdihs.rbcA[i] = old[i];
374                 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS+i];
375             }
376             break;
377         case F_CBTDIHS:
378             for (i = 0; (i < NR_CBTDIHS); i++)
379             {
380                 newparam->cbtdihs.cbtcA[i] = old[i];
381             }
382             break;
383         case F_FOURDIHS:
384             /* Read the dihedral parameters to temporary arrays,
385              * and convert them to the computationally faster
386              * Ryckaert-Bellemans form.
387              */
388             /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
389             newparam->rbdihs.rbcA[0] = old[1]+0.5*(old[0]+old[2]);
390             newparam->rbdihs.rbcA[1] = 0.5*(3.0*old[2]-old[0]);
391             newparam->rbdihs.rbcA[2] = 4.0*old[3]-old[1];
392             newparam->rbdihs.rbcA[3] = -2.0*old[2];
393             newparam->rbdihs.rbcA[4] = -4.0*old[3];
394             newparam->rbdihs.rbcA[5] = 0.0;
395
396             newparam->rbdihs.rbcB[0] = old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
397             newparam->rbdihs.rbcB[1] = 0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
398             newparam->rbdihs.rbcB[2] = 4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
399             newparam->rbdihs.rbcB[3] = -2.0*old[NR_FOURDIHS+2];
400             newparam->rbdihs.rbcB[4] = -4.0*old[NR_FOURDIHS+3];
401             newparam->rbdihs.rbcB[5] = 0.0;
402             break;
403         case F_CONSTR:
404         case F_CONSTRNC:
405             newparam->constr.dA = old[0];
406             newparam->constr.dB = old[1];
407             break;
408         case F_SETTLE:
409             newparam->settle.doh = old[0];
410             newparam->settle.dhh = old[1];
411             break;
412         case F_VSITE2:
413         case F_VSITE3:
414         case F_VSITE3FD:
415         case F_VSITE3OUT:
416         case F_VSITE4FD:
417         case F_VSITE4FDN:
418             newparam->vsite.a = old[0];
419             newparam->vsite.b = old[1];
420             newparam->vsite.c = old[2];
421             newparam->vsite.d = old[3];
422             newparam->vsite.e = old[4];
423             newparam->vsite.f = old[5];
424             break;
425         case F_VSITE3FAD:
426             newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
427             newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
428             newparam->vsite.c = old[2];
429             newparam->vsite.d = old[3];
430             newparam->vsite.e = old[4];
431             newparam->vsite.f = old[5];
432             break;
433         case F_VSITEN:
434             newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
435             newparam->vsiten.a = old[1];
436             break;
437         case F_CMAP:
438             newparam->cmap.cmapA = old[0];
439             newparam->cmap.cmapB = old[1];
440             break;
441         case F_GB12:
442         case F_GB13:
443         case F_GB14:
444             newparam->gb.sar  = old[0];
445             newparam->gb.st   = old[1];
446             newparam->gb.pi   = old[2];
447             newparam->gb.gbr  = old[3];
448             newparam->gb.bmlt = old[4];
449             break;
450         default:
451             gmx_fatal(FARGS, "unknown function type %d in %s line %d",
452                       ftype, __FILE__, __LINE__);
453     }
454     return 0;
455 }
456
457 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
458                         real forceparams[MAXFORCEPARAM], int comb, real reppow,
459                         int start, gmx_bool bAppend)
460 {
461     t_iparams newparam;
462     int       type;
463     int       rc;
464
465     if ( (rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
466     {
467         /* -1 means this interaction is all-zero and should not be added */
468         return rc;
469     }
470
471     if (!bAppend)
472     {
473         for (type = start; (type < ffparams->ntypes); type++)
474         {
475             if (ffparams->functype[type] == ftype)
476             {
477                 if (F_GB13 == ftype)
478                 {
479                     /* Occasionally, the way the 1-3 reference distance is
480                      * computed can lead to non-binary-identical results, but I
481                      * don't know why. */
482                     if ((gmx_within_tol(newparam.gb.sar,  ffparams->iparams[type].gb.sar,  1e-6)) &&
483                         (gmx_within_tol(newparam.gb.st,   ffparams->iparams[type].gb.st,   1e-6)) &&
484                         (gmx_within_tol(newparam.gb.pi,   ffparams->iparams[type].gb.pi,   1e-6)) &&
485                         (gmx_within_tol(newparam.gb.gbr,  ffparams->iparams[type].gb.gbr,  1e-6)) &&
486                         (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6)))
487                     {
488                         return type;
489                     }
490                 }
491                 else
492                 {
493                     if (memcmp(&newparam, &ffparams->iparams[type], (size_t)sizeof(newparam)) == 0)
494                     {
495                         return type;
496                     }
497                 }
498             }
499         }
500     }
501     else
502     {
503         type = ffparams->ntypes;
504     }
505     if (debug)
506     {
507         fprintf(debug, "copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
508                 type, ffparams->ntypes);
509     }
510     memcpy(&ffparams->iparams[type], &newparam, (size_t)sizeof(newparam));
511
512     ffparams->ntypes++;
513     ffparams->functype[type] = ftype;
514
515     return type;
516 }
517
518 static void append_interaction(t_ilist *ilist,
519                                int type, int nral, atom_id a[MAXATOMLIST])
520 {
521     int i, where1;
522
523     where1     = ilist->nr;
524     ilist->nr += nral+1;
525
526     ilist->iatoms[where1++] = type;
527     for (i = 0; (i < nral); i++)
528     {
529         ilist->iatoms[where1++] = a[i];
530     }
531 }
532
533 static void enter_function(t_params *p, t_functype ftype, int comb, real reppow,
534                            gmx_ffparams_t *ffparams, t_ilist *il,
535                            int *maxtypes,
536                            gmx_bool bNB, gmx_bool bAppend)
537 {
538     int     k, type, nr, nral, delta, start;
539
540     start = ffparams->ntypes;
541     nr    = p->nr;
542
543     for (k = 0; k < nr; k++)
544     {
545         if (*maxtypes <= ffparams->ntypes)
546         {
547             *maxtypes += 1000;
548             srenew(ffparams->functype, *maxtypes);
549             srenew(ffparams->iparams, *maxtypes);
550             if (debug)
551             {
552                 fprintf(debug, "%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
553                         __FILE__, __LINE__, *maxtypes);
554             }
555         }
556         type = enter_params(ffparams, ftype, p->param[k].c, comb, reppow, start, bAppend);
557         /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
558         if (!bNB && type >= 0)
559         {
560             nral  = NRAL(ftype);
561             delta = nr*(nral+1);
562             srenew(il->iatoms, il->nr+delta);
563             append_interaction(il, type, nral, p->param[k].a);
564         }
565     }
566 }
567
568 void convert_params(int atnr, t_params nbtypes[],
569                     t_molinfo *mi, int comb, double reppow, real fudgeQQ,
570                     gmx_mtop_t *mtop)
571 {
572     int             i, j, maxtypes, mt;
573     unsigned long   flags;
574     gmx_ffparams_t *ffp;
575     gmx_moltype_t  *molt;
576     t_params       *plist;
577
578     maxtypes = 0;
579
580     ffp           = &mtop->ffparams;
581     ffp->ntypes   = 0;
582     ffp->atnr     = atnr;
583     ffp->functype = NULL;
584     ffp->iparams  = NULL;
585     ffp->reppow   = reppow;
586
587     enter_function(&(nbtypes[F_LJ]),  (t_functype)F_LJ,    comb, reppow, ffp, NULL,
588                    &maxtypes, TRUE, TRUE);
589     enter_function(&(nbtypes[F_BHAM]), (t_functype)F_BHAM,  comb, reppow, ffp, NULL,
590                    &maxtypes, TRUE, TRUE);
591
592     for (mt = 0; mt < mtop->nmoltype; mt++)
593     {
594         molt = &mtop->moltype[mt];
595         for (i = 0; (i < F_NRE); i++)
596         {
597             molt->ilist[i].nr     = 0;
598             molt->ilist[i].iatoms = NULL;
599
600             plist = mi[mt].plist;
601
602             flags = interaction_function[i].flags;
603             if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
604                                                  (flags & IF_VSITE) ||
605                                                  (flags & IF_CONSTRAINT)))
606             {
607                 enter_function(&(plist[i]), (t_functype)i, comb, reppow,
608                                ffp, &molt->ilist[i],
609                                &maxtypes, FALSE, (i == F_POSRES  || i == F_FBPOSRES));
610             }
611         }
612     }
613     if (debug)
614     {
615         fprintf(debug, "%s, line %d: There are %d functypes in idef\n",
616                 __FILE__, __LINE__, ffp->ntypes);
617     }
618
619     ffp->fudgeQQ = fudgeQQ;
620 }