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