Project

General

Profile

forcerec.c

fix - Berk Hess, 10/12/2010 09:39 AM

 
1
/* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2
 *
3
 * 
4
 *                This source code is part of
5
 * 
6
 *                 G   R   O   M   A   C   S
7
 * 
8
 *          GROningen MAchine for Chemical Simulations
9
 * 
10
 *                        VERSION 3.2.0
11
 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12
 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13
 * Copyright (c) 2001-2004, The GROMACS development team,
14
 * check out http://www.gromacs.org for more information.
15

16
 * This program is free software; you can redistribute it and/or
17
 * modify it under the terms of the GNU General Public License
18
 * as published by the Free Software Foundation; either version 2
19
 * of the License, or (at your option) any later version.
20
 * 
21
 * If you want to redistribute modifications, please consider that
22
 * scientific software is very special. Version control is crucial -
23
 * bugs must be traceable. We will be happy to consider code for
24
 * inclusion in the official distribution, but derived work must not
25
 * be called official GROMACS. Details are found in the README & COPYING
26
 * files - if they are missing, get the official version at www.gromacs.org.
27
 * 
28
 * To help us fund GROMACS development, we humbly ask that you cite
29
 * the papers on the package - you can find them in the top README file.
30
 * 
31
 * For more info, check our website at http://www.gromacs.org
32
 * 
33
 * And Hey:
34
 * GROwing Monsters And Cloning Shrimps
35
 */
36
#ifdef HAVE_CONFIG_H
37
#include <config.h>
38
#endif
39

    
40
#include <math.h>
41
#include <string.h>
42
#include "sysstuff.h"
43
#include "typedefs.h"
44
#include "macros.h"
45
#include "smalloc.h"
46
#include "macros.h"
47
#include "physics.h"
48
#include "force.h"
49
#include "nonbonded.h"
50
#include "invblock.h"
51
#include "names.h"
52
#include "network.h"
53
#include "pbc.h"
54
#include "ns.h"
55
#include "mshift.h"
56
#include "txtdump.h"
57
#include "coulomb.h"
58
#include "mdrun.h"
59
#include "domdec.h"
60
#include "partdec.h"
61
#include "qmmm.h"
62
#include "copyrite.h"
63
#include "mtop_util.h"
64

    
65

    
66
#ifdef _MSC_VER
67
/* MSVC definition for __cpuid() */
68
#include <intrin.h>
69
#endif
70

    
71

    
72

    
73
t_forcerec *mk_forcerec(void)
74
{
75
  t_forcerec *fr;
76
  
77
  snew(fr,1);
78
  
79
  return fr;
80
}
81

    
82
#ifdef DEBUG
83
static void pr_nbfp(FILE *fp,real *nbfp,gmx_bool bBHAM,int atnr)
84
{
85
  int i,j;
86
  
87
  for(i=0; (i<atnr); i++) {
88
    for(j=0; (j<atnr); j++) {
89
      fprintf(fp,"%2d - %2d",i,j);
90
      if (bBHAM)
91
        fprintf(fp,"  a=%10g, b=%10g, c=%10g\n",BHAMA(nbfp,atnr,i,j),
92
                BHAMB(nbfp,atnr,i,j),BHAMC(nbfp,atnr,i,j));
93
      else
94
        fprintf(fp,"  c6=%10g, c12=%10g\n",C6(nbfp,atnr,i,j),
95
                C12(nbfp,atnr,i,j));
96
    }
97
  }
98
}
99
#endif
100

    
101
static real *mk_nbfp(const gmx_ffparams_t *idef,gmx_bool bBHAM)
102
{
103
  real *nbfp;
104
  int  i,j,k,atnr;
105
  
106
  atnr=idef->atnr;
107
  if (bBHAM) {
108
    snew(nbfp,3*atnr*atnr);
109
    for(i=k=0; (i<atnr); i++) {
110
      for(j=0; (j<atnr); j++,k++) {
111
        BHAMA(nbfp,atnr,i,j) = idef->iparams[k].bham.a;
112
        BHAMB(nbfp,atnr,i,j) = idef->iparams[k].bham.b;
113
        BHAMC(nbfp,atnr,i,j) = idef->iparams[k].bham.c;
114
      }
115
    }
116
  }
117
  else {
118
    snew(nbfp,2*atnr*atnr);
119
    for(i=k=0; (i<atnr); i++) {
120
      for(j=0; (j<atnr); j++,k++) {
121
        C6(nbfp,atnr,i,j)   = idef->iparams[k].lj.c6;
122
        C12(nbfp,atnr,i,j)  = idef->iparams[k].lj.c12;
123
      }
124
    }
125
  }
126
  return nbfp;
127
}
128

    
129
/* This routine sets fr->solvent_opt to the most common solvent in the 
130
 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in 
131
 * the fr->solvent_type array with the correct type (or esolNO).
132
 *
133
 * Charge groups that fulfill the conditions but are not identical to the
134
 * most common one will be marked as esolNO in the solvent_type array. 
135
 *
136
 * TIP3p is identical to SPC for these purposes, so we call it
137
 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
138
 * 
139
 * NOTE: QM particle should not
140
 * become an optimized solvent. Not even if there is only one charge
141
 * group in the Qm 
142
 */
143

    
144
typedef struct 
145
{
146
    int    model;          
147
    int    count;
148
    int    vdwtype[4];
149
    real   charge[4];
150
} solvent_parameters_t;
151

    
152
static void
153
check_solvent_cg(const gmx_moltype_t   *molt,
154
                 int                   cg0,
155
                 int                   nmol,
156
                 const unsigned char   *qm_grpnr,
157
                 const t_grps          *qm_grps,
158
                 t_forcerec *          fr,
159
                 int                   *n_solvent_parameters,
160
                 solvent_parameters_t  **solvent_parameters_p,
161
                 int                   cginfo,
162
                 int                   *cg_sp)
163
{
164
    const t_blocka *  excl;
165
    t_atom            *atom;
166
    int               j,k;
167
    int               j0,j1,nj;
168
    gmx_bool              perturbed;
169
    gmx_bool              has_vdw[4];
170
    gmx_bool              match;
171
    real              tmp_charge[4];
172
    int               tmp_vdwtype[4];
173
    int               tjA;
174
    gmx_bool              qm;
175
    solvent_parameters_t *solvent_parameters;
176

    
177
    /* We use a list with parameters for each solvent type. 
178
     * Every time we discover a new molecule that fulfills the basic 
179
     * conditions for a solvent we compare with the previous entries
180
     * in these lists. If the parameters are the same we just increment
181
     * the counter for that type, and otherwise we create a new type
182
     * based on the current molecule.
183
     *
184
     * Once we've finished going through all molecules we check which
185
     * solvent is most common, and mark all those molecules while we
186
     * clear the flag on all others.
187
     */   
188

    
189
    solvent_parameters = *solvent_parameters_p;
190

    
191
    /* Mark the cg first as non optimized */
192
    *cg_sp = -1;
193
    
194
    /* Check if this cg has no exclusions with atoms in other charge groups
195
     * and all atoms inside the charge group excluded.
196
     * We only have 3 or 4 atom solvent loops.
197
     */
198
    if (GET_CGINFO_EXCL_INTER(cginfo) ||
199
        !GET_CGINFO_EXCL_INTRA(cginfo))
200
    {
201
        return;
202
    }
203

    
204
    /* Get the indices of the first atom in this charge group */
205
    j0     = molt->cgs.index[cg0];
206
    j1     = molt->cgs.index[cg0+1];
207
    
208
    /* Number of atoms in our molecule */
209
    nj     = j1 - j0;
210

    
211
    if (debug) {
212
        fprintf(debug,
213
                "Moltype '%s': there are %d atoms in this charge group\n",
214
                *molt->name,nj);
215
    }
216
    
217
    /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
218
     * otherwise skip it.
219
     */
220
    if (nj<3 || nj>4)
221
    {
222
        return;
223
    }
224
    
225
    /* Check if we are doing QM on this group */
226
    qm = FALSE; 
227
    if (qm_grpnr != NULL)
228
    {
229
        for(j=j0 ; j<j1 && !qm; j++)
230
        {
231
            qm = (qm_grpnr[j] < qm_grps->nr - 1);
232
        }
233
    }
234
    /* Cannot use solvent optimization with QM */
235
    if (qm)
236
    {
237
        return;
238
    }
239
    
240
    atom = molt->atoms.atom;
241

    
242
    /* Still looks like a solvent, time to check parameters */
243
    
244
    /* If it is perturbed (free energy) we can't use the solvent loops,
245
     * so then we just skip to the next molecule.
246
     */   
247
    perturbed = FALSE; 
248
    
249
    for(j=j0; j<j1 && !perturbed; j++)
250
    {
251
        perturbed = PERTURBED(atom[j]);
252
    }
253
    
254
    if (perturbed)
255
    {
256
        return;
257
    }
258
    
259
    /* Now it's only a question if the VdW and charge parameters 
260
     * are OK. Before doing the check we compare and see if they are 
261
     * identical to a possible previous solvent type.
262
     * First we assign the current types and charges.    
263
     */
264
    for(j=0; j<nj; j++)
265
    {
266
        tmp_vdwtype[j] = atom[j0+j].type;
267
        tmp_charge[j]  = atom[j0+j].q;
268
    } 
269
    
270
    /* Does it match any previous solvent type? */
271
    for(k=0 ; k<*n_solvent_parameters; k++)
272
    {
273
        match = TRUE;
274
        
275
        
276
        /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
277
        if( (solvent_parameters[k].model==esolSPC   && nj!=3)  ||
278
            (solvent_parameters[k].model==esolTIP4P && nj!=4) )
279
            match = FALSE;
280
        
281
        /* Check that types & charges match for all atoms in molecule */
282
        for(j=0 ; j<nj && match==TRUE; j++)
283
        {                        
284
            if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
285
            {
286
                match = FALSE;
287
            }
288
            if(tmp_charge[j] != solvent_parameters[k].charge[j])
289
            {
290
                match = FALSE;
291
            }
292
        }
293
        if (match == TRUE)
294
        {
295
            /* Congratulations! We have a matched solvent.
296
             * Flag it with this type for later processing.
297
             */
298
            *cg_sp = k;
299
            solvent_parameters[k].count += nmol;
300

    
301
            /* We are done with this charge group */
302
            return;
303
        }
304
    }
305
    
306
    /* If we get here, we have a tentative new solvent type.
307
     * Before we add it we must check that it fulfills the requirements
308
     * of the solvent optimized loops. First determine which atoms have
309
     * VdW interactions.   
310
     */
311
    for(j=0; j<nj; j++) 
312
    {
313
        has_vdw[j] = FALSE;
314
        tjA        = tmp_vdwtype[j];
315
        
316
        /* Go through all other tpes and see if any have non-zero
317
         * VdW parameters when combined with this one.
318
         */   
319
        for(k=0; k<fr->ntype && (has_vdw[j]==FALSE); k++)
320
        {
321
            /* We already checked that the atoms weren't perturbed,
322
             * so we only need to check state A now.
323
             */ 
324
            if (fr->bBHAM) 
325
            {
326
                has_vdw[j] = (has_vdw[j] || 
327
                              (BHAMA(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
328
                              (BHAMB(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
329
                              (BHAMC(fr->nbfp,fr->ntype,tjA,k) != 0.0));
330
            }
331
            else
332
            {
333
                /* Standard LJ */
334
                has_vdw[j] = (has_vdw[j] || 
335
                              (C6(fr->nbfp,fr->ntype,tjA,k)  != 0.0) ||
336
                              (C12(fr->nbfp,fr->ntype,tjA,k) != 0.0));
337
            }
338
        }
339
    }
340
    
341
    /* Now we know all we need to make the final check and assignment. */
342
    if (nj == 3)
343
    {
344
        /* So, is it an SPC?
345
         * For this we require thatn all atoms have charge, 
346
         * the charges on atom 2 & 3 should be the same, and only
347
         * atom 1 should have VdW.
348
         */
349
        if (has_vdw[0] == TRUE && 
350
            has_vdw[1] == FALSE &&
351
            has_vdw[2] == FALSE &&
352
            tmp_charge[0]  != 0 &&
353
            tmp_charge[1]  != 0 &&
354
            tmp_charge[2]  == tmp_charge[1])
355
        {
356
            srenew(solvent_parameters,*n_solvent_parameters+1);
357
            solvent_parameters[*n_solvent_parameters].model = esolSPC;
358
            solvent_parameters[*n_solvent_parameters].count = nmol;
359
            for(k=0;k<3;k++)
360
            {
361
                solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
362
                solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
363
            }
364

    
365
            *cg_sp = *n_solvent_parameters;
366
            (*n_solvent_parameters)++;
367
        }
368
    }
369
    else if (nj==4)
370
    {
371
        /* Or could it be a TIP4P?
372
         * For this we require thatn atoms 2,3,4 have charge, but not atom 1. 
373
         * Only atom 1 should have VdW.
374
         */
375
        if(has_vdw[0] == TRUE && 
376
           has_vdw[1] == FALSE &&
377
           has_vdw[2] == FALSE &&
378
           has_vdw[3] == FALSE &&
379
           tmp_charge[0]  == 0 &&
380
           tmp_charge[1]  != 0 &&
381
           tmp_charge[2]  == tmp_charge[1] &&
382
           tmp_charge[3]  != 0)
383
        {
384
            srenew(solvent_parameters,*n_solvent_parameters+1);
385
            solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
386
            solvent_parameters[*n_solvent_parameters].count = nmol;
387
            for(k=0;k<4;k++)
388
            {
389
                solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
390
                solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
391
            }
392
            
393
            *cg_sp = *n_solvent_parameters;
394
            (*n_solvent_parameters)++;
395
        }
396
    }
397

    
398
    *solvent_parameters_p = solvent_parameters;
399
}
400

    
401
static void
402
check_solvent(FILE *                fp,
403
              const gmx_mtop_t *    mtop,
404
              t_forcerec *          fr,
405
              cginfo_mb_t           *cginfo_mb)
406
{
407
    const t_block *   cgs;
408
    const t_block *   mols;
409
    const gmx_moltype_t *molt;
410
    int               mb,mol,cg_mol,at_offset,cg_offset,am,cgm,i,nmol_ch,nmol;
411
    int               n_solvent_parameters;
412
    solvent_parameters_t *solvent_parameters;
413
    int               **cg_sp;
414
    int               bestsp,bestsol;
415

    
416
    if (debug)
417
    {
418
        fprintf(debug,"Going to determine what solvent types we have.\n");
419
    }
420

    
421
    mols = &mtop->mols;
422

    
423
    n_solvent_parameters = 0;
424
    solvent_parameters = NULL;
425
    /* Allocate temporary array for solvent type */
426
    snew(cg_sp,mtop->nmolblock);
427

    
428
    cg_offset = 0;
429
    at_offset = 0;
430
    for(mb=0; mb<mtop->nmolblock; mb++)
431
    {
432
        molt = &mtop->moltype[mtop->molblock[mb].type];
433
        cgs  = &molt->cgs;
434
        /* Here we have to loop over all individual molecules
435
         * because we need to check for QMMM particles.
436
         */
437
        snew(cg_sp[mb],cginfo_mb[mb].cg_mod);
438
        nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
439
        nmol    = mtop->molblock[mb].nmol/nmol_ch;
440
        for(mol=0; mol<nmol_ch; mol++)
441
        {
442
            cgm = mol*cgs->nr;
443
            am  = mol*cgs->index[cgs->nr];
444
            for(cg_mol=0; cg_mol<cgs->nr; cg_mol++)
445
            {
446
                check_solvent_cg(molt,cg_mol,nmol,
447
                                 mtop->groups.grpnr[egcQMMM] ?
448
                                 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
449
                                 &mtop->groups.grps[egcQMMM],
450
                                 fr,
451
                                 &n_solvent_parameters,&solvent_parameters,
452
                                 cginfo_mb[mb].cginfo[cgm+cg_mol],
453
                                 &cg_sp[mb][cgm+cg_mol]);
454
            }
455
        }
456
        cg_offset += cgs->nr;
457
        at_offset += cgs->index[cgs->nr];
458
    }
459

    
460
    /* Puh! We finished going through all charge groups.
461
     * Now find the most common solvent model.
462
     */   
463
    
464
    /* Most common solvent this far */
465
    bestsp = -2;
466
    for(i=0;i<n_solvent_parameters;i++)
467
    {
468
        if (bestsp == -2 ||
469
            solvent_parameters[i].count > solvent_parameters[bestsp].count)
470
        {
471
            bestsp = i;
472
        }
473
    }
474
    
475
    if (bestsp >= 0)
476
    {
477
        bestsol = solvent_parameters[bestsp].model;
478
    }
479
    else
480
    {
481
        bestsol = esolNO;
482
    }
483
    
484
#ifdef DISABLE_WATER_NLIST
485
        bestsol = esolNO;
486
#endif
487

    
488
    fr->nWatMol = 0;
489
    for(mb=0; mb<mtop->nmolblock; mb++)
490
    {
491
        cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
492
        nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
493
        for(i=0; i<cginfo_mb[mb].cg_mod; i++)
494
        {
495
            if (cg_sp[mb][i] == bestsp)
496
            {
497
                SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],bestsol);
498
                fr->nWatMol += nmol;
499
            }
500
            else
501
            {
502
                SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],esolNO);
503
            }
504
        }
505
        sfree(cg_sp[mb]);
506
    }
507
    sfree(cg_sp);
508
    
509
    if (bestsol != esolNO && fp!=NULL)
510
    {
511
        fprintf(fp,"\nEnabling %s-like water optimization for %d molecules.\n\n",
512
                esol_names[bestsol],
513
                solvent_parameters[bestsp].count);
514
    }
515

    
516
    sfree(solvent_parameters);
517
    fr->solvent_opt = bestsol;
518
}
519

    
520
static cginfo_mb_t *init_cginfo_mb(FILE *fplog,const gmx_mtop_t *mtop,
521
                                   t_forcerec *fr,gmx_bool bNoSolvOpt,
522
                                   gmx_bool *bExcl_IntraCGAll_InterCGNone)
523
{
524
    const t_block *cgs;
525
    const t_blocka *excl;
526
    const gmx_moltype_t *molt;
527
    const gmx_molblock_t *molb;
528
    cginfo_mb_t *cginfo_mb;
529
    int  *cginfo;
530
    int  cg_offset,a_offset,cgm,am;
531
    int  mb,m,ncg_tot,cg,a0,a1,gid,ai,j,aj,excl_nalloc;
532
    gmx_bool bId,*bExcl,bExclIntraAll,bExclInter;
533

    
534
    ncg_tot = ncg_mtop(mtop);
535
    snew(cginfo_mb,mtop->nmolblock);
536

    
537
    *bExcl_IntraCGAll_InterCGNone = TRUE;
538

    
539
    excl_nalloc = 10;
540
    snew(bExcl,excl_nalloc);
541
    cg_offset = 0;
542
    a_offset  = 0;
543
    for(mb=0; mb<mtop->nmolblock; mb++)
544
    {
545
        molb = &mtop->molblock[mb];
546
        molt = &mtop->moltype[molb->type];
547
        cgs  = &molt->cgs;
548
        excl = &molt->excls;
549

    
550
        /* Check if the cginfo is identical for all molecules in this block.
551
         * If so, we only need an array of the size of one molecule.
552
         * Otherwise we make an array of #mol times #cgs per molecule.
553
         */
554
        bId = TRUE;
555
        am = 0;
556
        for(m=0; m<molb->nmol; m++)
557
        {
558
            am = m*cgs->index[cgs->nr];
559
            for(cg=0; cg<cgs->nr; cg++)
560
            {
561
                a0 = cgs->index[cg];
562
                a1 = cgs->index[cg+1];
563
                if (ggrpnr(&mtop->groups,egcENER,a_offset+am+a0) !=
564
                    ggrpnr(&mtop->groups,egcENER,a_offset   +a0))
565
                {
566
                    bId = FALSE;
567
                }
568
                if (mtop->groups.grpnr[egcQMMM] != NULL)
569
                {
570
                    for(ai=a0; ai<a1; ai++)
571
                    {
572
                        if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
573
                            mtop->groups.grpnr[egcQMMM][a_offset   +ai])
574
                        {
575
                            bId = FALSE;
576
                        }
577
                    }
578
                }
579
            }
580
        }
581

    
582
        cginfo_mb[mb].cg_start = cg_offset;
583
        cginfo_mb[mb].cg_end   = cg_offset + molb->nmol*cgs->nr;
584
        cginfo_mb[mb].cg_mod   = (bId ? 1 : molb->nmol)*cgs->nr;
585
        snew(cginfo_mb[mb].cginfo,cginfo_mb[mb].cg_mod);
586
        cginfo = cginfo_mb[mb].cginfo;
587

    
588
        for(m=0; m<(bId ? 1 : molb->nmol); m++)
589
        {
590
            cgm = m*cgs->nr;
591
            am  = m*cgs->index[cgs->nr];
592
            for(cg=0; cg<cgs->nr; cg++)
593
            {
594
                a0 = cgs->index[cg];
595
                a1 = cgs->index[cg+1];
596

    
597
                /* Store the energy group in cginfo */
598
                gid = ggrpnr(&mtop->groups,egcENER,a_offset+am+a0);
599
                SET_CGINFO_GID(cginfo[cgm+cg],gid);
600
                
601
                /* Check the intra/inter charge group exclusions */
602
                if (a1-a0 > excl_nalloc) {
603
                    excl_nalloc = a1 - a0;
604
                    srenew(bExcl,excl_nalloc);
605
                }
606
                /* bExclIntraAll: all intra cg interactions excluded
607
                 * bExclInter:    any inter cg interactions excluded
608
                 */
609
                bExclIntraAll = TRUE;
610
                bExclInter    = FALSE;
611
                for(ai=a0; ai<a1; ai++) {
612
                    /* Clear the exclusion list for atom ai */
613
                    for(aj=a0; aj<a1; aj++) {
614
                        bExcl[aj-a0] = FALSE;
615
                    }
616
                    /* Loop over all the exclusions of atom ai */
617
                    for(j=excl->index[ai]; j<excl->index[ai+1]; j++)
618
                    {
619
                        aj = excl->a[j];
620
                        if (aj < a0 || aj >= a1)
621
                        {
622
                            bExclInter = TRUE;
623
                        }
624
                        else
625
                        {
626
                            bExcl[aj-a0] = TRUE;
627
                        }
628
                    }
629
                    /* Check if ai excludes a0 to a1 */
630
                    for(aj=a0; aj<a1; aj++)
631
                    {
632
                        if (!bExcl[aj-a0])
633
                        {
634
                            bExclIntraAll = FALSE;
635
                        }
636
                    }
637
                }
638
                if (bExclIntraAll)
639
                {
640
                    SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
641
                }
642
                if (bExclInter)
643
                {
644
                    SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
645
                }
646
                if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
647
                {
648
                    /* The size in cginfo is currently only read with DD */
649
                    gmx_fatal(FARGS,"A charge group has size %d which is larger than the limit of %d atoms",a1-a0,MAX_CHARGEGROUP_SIZE);
650
                }
651
                SET_CGINFO_NATOMS(cginfo[cgm+cg],a1-a0);
652

    
653
                if (!bExclIntraAll || bExclInter)
654
                {
655
                    *bExcl_IntraCGAll_InterCGNone = FALSE;
656
                }
657
            }
658
        }
659
        cg_offset += molb->nmol*cgs->nr;
660
        a_offset  += molb->nmol*cgs->index[cgs->nr];
661
    }
662
    sfree(bExcl);
663
    
664
    /* the solvent optimizer is called after the QM is initialized,
665
     * because we don't want to have the QM subsystemto become an
666
     * optimized solvent
667
     */
668

    
669
    check_solvent(fplog,mtop,fr,cginfo_mb);
670
    
671
    if (getenv("GMX_NO_SOLV_OPT"))
672
    {
673
        if (fplog)
674
        {
675
            fprintf(fplog,"Found environment variable GMX_NO_SOLV_OPT.\n"
676
                    "Disabling all solvent optimization\n");
677
        }
678
        fr->solvent_opt = esolNO;
679
    }
680
    if (bNoSolvOpt)
681
    {
682
        fr->solvent_opt = esolNO;
683
    }
684
    if (!fr->solvent_opt)
685
    {
686
        for(mb=0; mb<mtop->nmolblock; mb++)
687
        {
688
            for(cg=0; cg<cginfo_mb[mb].cg_mod; cg++)
689
            {
690
                SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg],esolNO);
691
            }
692
        }
693
    }
694
    
695
    return cginfo_mb;
696
}
697

    
698
static int *cginfo_expand(int nmb,cginfo_mb_t *cgi_mb)
699
{
700
    int ncg,mb,cg;
701
    int *cginfo;
702

    
703
    ncg = cgi_mb[nmb-1].cg_end;
704
    snew(cginfo,ncg);
705
    mb = 0;
706
    for(cg=0; cg<ncg; cg++)
707
    {
708
        while (cg >= cgi_mb[mb].cg_end)
709
        {
710
            mb++;
711
        }
712
        cginfo[cg] =
713
            cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
714
    }
715

    
716
    return cginfo;
717
}
718

    
719
static void set_chargesum(FILE *log,t_forcerec *fr,const gmx_mtop_t *mtop)
720
{
721
    double qsum;
722
    int    mb,nmol,i;
723
    const t_atoms *atoms;
724
    
725
    qsum = 0;
726
    for(mb=0; mb<mtop->nmolblock; mb++)
727
    {
728
        nmol  = mtop->molblock[mb].nmol;
729
        atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
730
        for(i=0; i<atoms->nr; i++)
731
        {
732
            qsum += nmol*atoms->atom[i].q;
733
        }
734
    }
735
    fr->qsum[0] = qsum;
736
    if (fr->efep != efepNO)
737
    {
738
        qsum = 0;
739
        for(mb=0; mb<mtop->nmolblock; mb++)
740
        {
741
            nmol  = mtop->molblock[mb].nmol;
742
            atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
743
            for(i=0; i<atoms->nr; i++)
744
            {
745
                qsum += nmol*atoms->atom[i].qB;
746
            }
747
            fr->qsum[1] = qsum;
748
        }
749
    }
750
    else
751
    {
752
        fr->qsum[1] = fr->qsum[0];
753
    }
754
    if (log) {
755
        if (fr->efep == efepNO)
756
            fprintf(log,"System total charge: %.3f\n",fr->qsum[0]);
757
        else
758
            fprintf(log,"System total charge, top. A: %.3f top. B: %.3f\n",
759
                    fr->qsum[0],fr->qsum[1]);
760
    }
761
}
762

    
763
void update_forcerec(FILE *log,t_forcerec *fr,matrix box)
764
{
765
    if (fr->eeltype == eelGRF)
766
    {
767
        calc_rffac(NULL,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
768
                   fr->rcoulomb,fr->temp,fr->zsquare,box,
769
                   &fr->kappa,&fr->k_rf,&fr->c_rf);
770
    }
771
}
772

    
773
void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,const gmx_mtop_t *mtop)
774
{
775
    const t_atoms *atoms,*atoms_tpi;
776
    const t_blocka *excl;
777
    int    mb,nmol,nmolc,i,j,tpi,tpj,j1,j2,k,n,nexcl,q;
778
#if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)    
779
    long long int  npair,npair_ij,tmpi,tmpj;
780
#else
781
    double npair, npair_ij,tmpi,tmpj;
782
#endif
783
    double csix,ctwelve;
784
    int    ntp,*typecount;
785
    gmx_bool   bBHAM;
786
    real   *nbfp;
787

    
788
    ntp = fr->ntype;
789
    bBHAM = fr->bBHAM;
790
    nbfp = fr->nbfp;
791
    
792
    for(q=0; q<(fr->efep==efepNO ? 1 : 2); q++) {
793
        csix = 0;
794
        ctwelve = 0;
795
        npair = 0;
796
        nexcl = 0;
797
        if (!fr->n_tpi) {
798
            /* Count the types so we avoid natoms^2 operations */
799
            snew(typecount,ntp);
800
            for(mb=0; mb<mtop->nmolblock; mb++) {
801
                nmol  = mtop->molblock[mb].nmol;
802
                atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
803
                for(i=0; i<atoms->nr; i++) {
804
                    if (q == 0)
805
                    {
806
                        tpi = atoms->atom[i].type;
807
                    }
808
                    else
809
                    {
810
                        tpi = atoms->atom[i].typeB;
811
                    }
812
                    typecount[tpi] += nmol;
813
                }
814
            }
815
            for(tpi=0; tpi<ntp; tpi++) {
816
                for(tpj=tpi; tpj<ntp; tpj++) {
817
                    tmpi = typecount[tpi];
818
                    tmpj = typecount[tpj];
819
                    if (tpi != tpj)
820
                    {
821
                        npair_ij = tmpi*tmpj;
822
                    }
823
                    else
824
                    {
825
                        npair_ij = tmpi*(tmpi - 1)/2;
826
                    }
827
                    if (bBHAM) {
828
                        csix    += npair_ij*BHAMC(nbfp,ntp,tpi,tpj);
829
                    } else {
830
                        csix    += npair_ij*   C6(nbfp,ntp,tpi,tpj);
831
                        ctwelve += npair_ij*  C12(nbfp,ntp,tpi,tpj);
832
                    }
833
                    npair += npair_ij;
834
                }
835
            }
836
            sfree(typecount);
837
            /* Subtract the excluded pairs.
838
             * The main reason for substracting exclusions is that in some cases
839
             * some combinations might never occur and the parameters could have
840
             * any value. These unused values should not influence the dispersion
841
             * correction.
842
             */
843
            for(mb=0; mb<mtop->nmolblock; mb++) {
844
                nmol  = mtop->molblock[mb].nmol;
845
                atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
846
                excl  = &mtop->moltype[mtop->molblock[mb].type].excls;
847
                for(i=0; (i<atoms->nr); i++) {
848
                    if (q == 0)
849
                    {
850
                        tpi = atoms->atom[i].type;
851
                    }
852
                    else
853
                    {
854
                        tpi = atoms->atom[i].typeB;
855
                    }
856
                    j1  = excl->index[i];
857
                    j2  = excl->index[i+1];
858
                    for(j=j1; j<j2; j++) {
859
                        k = excl->a[j];
860
                        if (k > i)
861
                        {
862
                            if (q == 0)
863
                            {
864
                                tpj = atoms->atom[k].type;
865
                            }
866
                            else
867
                            {
868
                                tpj = atoms->atom[k].typeB;
869
                            }
870
                            if (bBHAM) {
871
                                csix -= nmol*BHAMC(nbfp,ntp,tpi,tpj);
872
                            } else {
873
                                csix    -= nmol*C6 (nbfp,ntp,tpi,tpj);
874
                                ctwelve -= nmol*C12(nbfp,ntp,tpi,tpj);
875
                            }
876
                            nexcl += nmol;
877
                        }
878
                    }
879
                }
880
            }
881
        } else {
882
            /* Only correct for the interaction of the test particle
883
             * with the rest of the system.
884
             */
885
            atoms_tpi =
886
                &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
887

    
888
            npair = 0;
889
            for(mb=0; mb<mtop->nmolblock; mb++) {
890
                nmol  = mtop->molblock[mb].nmol;
891
                atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
892
                for(j=0; j<atoms->nr; j++) {
893
                    nmolc = nmol;
894
                    /* Remove the interaction of the test charge group
895
                     * with itself.
896
                     */
897
                    if (mb == mtop->nmolblock-1)
898
                    {
899
                        nmolc--;
900
                        
901
                        if (mb == 0 && nmol == 1)
902
                        {
903
                            gmx_fatal(FARGS,"Old format tpr with TPI, please generate a new tpr file");
904
                        }
905
                    }
906
                    if (q == 0)
907
                    {
908
                        tpj = atoms->atom[j].type;
909
                    }
910
                    else
911
                    {
912
                        tpj = atoms->atom[j].typeB;
913
                    }
914
                    for(i=0; i<fr->n_tpi; i++)
915
                    {
916
                        if (q == 0)
917
                        {
918
                            tpi = atoms_tpi->atom[i].type;
919
                        }
920
                        else
921
                        {
922
                            tpi = atoms_tpi->atom[i].typeB;
923
                        }
924
                        if (bBHAM)
925
                        {
926
                            csix    += nmolc*BHAMC(nbfp,ntp,tpi,tpj);
927
                        }
928
                        else
929
                        {
930
                            csix    += nmolc*C6 (nbfp,ntp,tpi,tpj);
931
                            ctwelve += nmolc*C12(nbfp,ntp,tpi,tpj);
932
                        }
933
                        npair += nmolc;
934
                    }
935
                }
936
            }
937
        }
938
        if (npair - nexcl <= 0 && fplog) {
939
            fprintf(fplog,"\nWARNING: There are no atom pairs for dispersion correction\n\n");
940
            csix     = 0;
941
            ctwelve  = 0;
942
        } else {
943
            csix    /= npair - nexcl;
944
            ctwelve /= npair - nexcl;
945
        }
946
        if (debug) {
947
            fprintf(debug,"Counted %d exclusions\n",nexcl);
948
            fprintf(debug,"Average C6 parameter is: %10g\n",(double)csix);
949
            fprintf(debug,"Average C12 parameter is: %10g\n",(double)ctwelve);
950
        }
951
        fr->avcsix[q]    = csix;
952
        fr->avctwelve[q] = ctwelve;
953
    }
954
    if (fplog != NULL)
955
    {
956
        if (fr->eDispCorr == edispcAllEner ||
957
            fr->eDispCorr == edispcAllEnerPres)
958
        {
959
            fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
960
                    fr->avcsix[0],fr->avctwelve[0]);
961
        }
962
        else
963
        {
964
            fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e\n",fr->avcsix[0]);
965
        }
966
    }
967
}
968

    
969

    
970
static void set_bham_b_max(FILE *fplog,t_forcerec *fr,
971
                           const gmx_mtop_t *mtop)
972
{
973
    const t_atoms *at1,*at2;
974
    int  mt1,mt2,i,j,tpi,tpj,ntypes;
975
    real b,bmin;
976
    real *nbfp;
977

    
978
    if (fplog)
979
    {
980
        fprintf(fplog,"Determining largest Buckingham b parameter for table\n");
981
    }
982
    nbfp   = fr->nbfp;
983
    ntypes = fr->ntype;
984
    
985
    bmin           = -1;
986
    fr->bham_b_max = 0;
987
    for(mt1=0; mt1<mtop->nmoltype; mt1++)
988
    {
989
        at1 = &mtop->moltype[mt1].atoms;
990
        for(i=0; (i<at1->nr); i++)
991
        {
992
            tpi = at1->atom[i].type;
993
            if (tpi >= ntypes)
994
                gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",i,tpi,ntypes);
995
            
996
            for(mt2=mt1; mt2<mtop->nmoltype; mt2++)
997
            {
998
                at2 = &mtop->moltype[mt2].atoms;
999
                for(j=0; (j<at2->nr); j++) {
1000
                    tpj = at2->atom[j].type;
1001
                    if (tpj >= ntypes)
1002
                    {
1003
                        gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",j,tpj,ntypes);
1004
                    }
1005
                    b = BHAMB(nbfp,ntypes,tpi,tpj);
1006
                    if (b > fr->bham_b_max)
1007
                    {
1008
                        fr->bham_b_max = b;
1009
                    }
1010
                    if ((b < bmin) || (bmin==-1))
1011
                    {
1012
                        bmin = b;
1013
                    }
1014
                }
1015
            }
1016
        }
1017
    }
1018
    if (fplog)
1019
    {
1020
        fprintf(fplog,"Buckingham b parameters, min: %g, max: %g\n",
1021
                bmin,fr->bham_b_max);
1022
    }
1023
}
1024

    
1025
static void make_nbf_tables(FILE *fp,const output_env_t oenv,
1026
                            t_forcerec *fr,real rtab,
1027
                            const t_commrec *cr,
1028
                            const char *tabfn,char *eg1,char *eg2,
1029
                            t_nblists *nbl)
1030
{
1031
  char buf[STRLEN];
1032
  int i,j;
1033
  void *      p_tmp;
1034

    
1035
  if (tabfn == NULL) {
1036
    if (debug)
1037
      fprintf(debug,"No table file name passed, can not read table, can not do non-bonded interactions\n");
1038
    return;
1039
  }
1040
    
1041
  sprintf(buf,"%s",tabfn);
1042
  if (eg1 && eg2)
1043
    /* Append the two energy group names */
1044
    sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,"_%s_%s.%s",
1045
            eg1,eg2,ftp2ext(efXVG));
1046
  nbl->tab = make_tables(fp,oenv,fr,MASTER(cr),buf,rtab,0);
1047
  /* Copy the contents of the table to separate coulomb and LJ tables too,
1048
   * to improve cache performance.
1049
   */
1050

    
1051
  /* For performance reasons we want
1052
   * the table data to be aligned to 16-byte. This is accomplished
1053
   * by allocating 16 bytes extra to a temporary pointer, and then
1054
   * calculating an aligned pointer. This new pointer must not be
1055
   * used in a free() call, but thankfully we're sloppy enough not
1056
   * to do this...
1057
   */
1058
  
1059
  /* 8 fp entries per vdw table point, n+1 points, and 16 bytes extra to align it. */
1060
  p_tmp = malloc(8*(nbl->tab.n+1)*sizeof(real)+16);
1061
  
1062
  /* align it - size_t has the same same as a pointer */
1063
  nbl->vdwtab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));  
1064

    
1065
  /* 4 fp entries per coul table point, n+1 points, and 16 bytes extra to align it. */
1066
  p_tmp = malloc(4*(nbl->tab.n+1)*sizeof(real)+16);
1067
  
1068
  /* align it - size_t has the same same as a pointer */
1069
  nbl->coultab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));  
1070

    
1071
  
1072
  for(i=0; i<=nbl->tab.n; i++) {
1073
    for(j=0; j<4; j++)
1074
      nbl->coultab[4*i+j] = nbl->tab.tab[12*i+j];
1075
    for(j=0; j<8; j++)
1076
      nbl->vdwtab [8*i+j] = nbl->tab.tab[12*i+4+j];
1077
  }
1078
}
1079

    
1080
static void count_tables(int ftype1,int ftype2,const gmx_mtop_t *mtop,
1081
                         int *ncount,int **count)
1082
{
1083
    const gmx_moltype_t *molt;
1084
    const t_ilist *il;
1085
    int mt,ftype,stride,i,j,tabnr;
1086
    
1087
    for(mt=0; mt<mtop->nmoltype; mt++)
1088
    {
1089
        molt = &mtop->moltype[mt];
1090
        for(ftype=0; ftype<F_NRE; ftype++)
1091
        {
1092
            if (ftype == ftype1 || ftype == ftype2) {
1093
                il = &molt->ilist[ftype];
1094
                stride = 1 + NRAL(ftype);
1095
                for(i=0; i<il->nr; i+=stride) {
1096
                    tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1097
                    if (tabnr < 0)
1098
                        gmx_fatal(FARGS,"A bonded table number is smaller than 0: %d\n",tabnr);
1099
                    if (tabnr >= *ncount) {
1100
                        srenew(*count,tabnr+1);
1101
                        for(j=*ncount; j<tabnr+1; j++)
1102
                            (*count)[j] = 0;
1103
                        *ncount = tabnr+1;
1104
                    }
1105
                    (*count)[tabnr]++;
1106
                }
1107
            }
1108
        }
1109
    }
1110
}
1111

    
1112
static bondedtable_t *make_bonded_tables(FILE *fplog,
1113
                                         int ftype1,int ftype2,
1114
                                         const gmx_mtop_t *mtop,
1115
                                         const char *basefn,const char *tabext)
1116
{
1117
    int  i,ncount,*count;
1118
    char tabfn[STRLEN];
1119
    bondedtable_t *tab;
1120
    
1121
    tab = NULL;
1122
    
1123
    ncount = 0;
1124
    count = NULL;
1125
    count_tables(ftype1,ftype2,mtop,&ncount,&count);
1126
    
1127
    if (ncount > 0) {
1128
        snew(tab,ncount);
1129
        for(i=0; i<ncount; i++) {
1130
            if (count[i] > 0) {
1131
                sprintf(tabfn,"%s",basefn);
1132
                sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1,"_%s%d.%s",
1133
                        tabext,i,ftp2ext(efXVG));
1134
                tab[i] = make_bonded_table(fplog,tabfn,NRAL(ftype1)-2);
1135
            }
1136
        }
1137
        sfree(count);
1138
    }
1139
  
1140
    return tab;
1141
}
1142

    
1143
void forcerec_set_ranges(t_forcerec *fr,
1144
                         int ncg_home,int ncg_force,
1145
                         int natoms_force,
1146
                         int natoms_force_constr,int natoms_f_novirsum)
1147
{
1148
    fr->cg0 = 0;
1149
    fr->hcg = ncg_home;
1150

    
1151
    /* fr->ncg_force is unused in the standard code,
1152
     * but it can be useful for modified code dealing with charge groups.
1153
     */
1154
    fr->ncg_force           = ncg_force;
1155
    fr->natoms_force        = natoms_force;
1156
    fr->natoms_force_constr = natoms_force_constr;
1157

    
1158
    if (fr->natoms_force_constr > fr->nalloc_force)
1159
    {
1160
        fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1161

    
1162
        if (fr->bTwinRange)
1163
        {
1164
            srenew(fr->f_twin,fr->nalloc_force);
1165
        }
1166
    }
1167

    
1168
    if (fr->bF_NoVirSum)
1169
    {
1170
        fr->f_novirsum_n = natoms_f_novirsum;
1171
        if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1172
        {
1173
            fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1174
            srenew(fr->f_novirsum_alloc,fr->f_novirsum_nalloc);
1175
        }
1176
    }
1177
    else
1178
    {
1179
        fr->f_novirsum_n = 0;
1180
    }
1181
}
1182

    
1183
static real cutoff_inf(real cutoff)
1184
{
1185
    if (cutoff == 0)
1186
    {
1187
        cutoff = GMX_CUTOFF_INF;
1188
    }
1189

    
1190
    return cutoff;
1191
}
1192

    
1193
gmx_bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
1194
                      gmx_bool bPrintNote,t_commrec *cr,FILE *fp)
1195
{
1196
    gmx_bool bAllvsAll;
1197

    
1198
    bAllvsAll =
1199
        (
1200
         ir->rlist==0            &&
1201
         ir->rcoulomb==0         &&
1202
         ir->rvdw==0             &&
1203
         ir->ePBC==epbcNONE      &&
1204
         ir->vdwtype==evdwCUT    &&
1205
         ir->coulombtype==eelCUT &&
1206
         ir->efep==efepNO        &&
1207
         (ir->implicit_solvent == eisNO || 
1208
          (ir->implicit_solvent==eisGBSA && (ir->gb_algorithm==egbSTILL || 
1209
                                             ir->gb_algorithm==egbHCT   || 
1210
                                             ir->gb_algorithm==egbOBC))) &&
1211
         getenv("GMX_NO_ALLVSALL") == NULL
1212
            );
1213
    
1214
    if (bAllvsAll && ir->opts.ngener > 1)
1215
    {
1216
        const char *note="NOTE: Can not use all-vs-all force loops, because there are multiple energy monitor groups; you might get significantly higher performance when using only a single energy monitor group.\n";
1217

    
1218
        if (bPrintNote)
1219
        {
1220
            if (MASTER(cr))
1221
            {
1222
                fprintf(stderr,"\n%s\n",note);
1223
            }
1224
            if (fp != NULL)
1225
            {
1226
                fprintf(fp,"\n%s\n",note);
1227
            }
1228
        }
1229
        bAllvsAll = FALSE;
1230
    }
1231

    
1232
    if(bAllvsAll && fp && MASTER(cr))
1233
    {
1234
        fprintf(fp,"\nUsing accelerated all-vs-all kernels.\n\n");
1235
    }
1236
    
1237
    return bAllvsAll;
1238
}
1239

    
1240

    
1241
/* Return 1 if SSE2 support is present, otherwise 0. */
1242
static int 
1243
forcerec_check_sse2()
1244
{
1245
#if ( defined(GMX_IA32_SSE2) || defined(GMX_X86_64_SSE2) || defined(GMX_IA32_SSE) || defined(GMX_X86_64_SSE)|| defined(GMX_SSE2) )
1246
        unsigned int level;
1247
        unsigned int _eax,_ebx,_ecx,_edx;
1248
        int status;
1249
        int CPUInfo[4];
1250
        
1251
        level = 1;
1252
#ifdef _MSC_VER
1253
        __cpuid(CPUInfo,1);
1254
        
1255
        _eax=CPUInfo[0];
1256
        _ebx=CPUInfo[1];
1257
        _ecx=CPUInfo[2];
1258
        _edx=CPUInfo[3];
1259
        
1260
#elif defined(__x86_64__)
1261
        /* GCC 64-bit inline asm */
1262
        __asm__ ("push %%rbx\n\tcpuid\n\tpop %%rbx\n"                 \
1263
                         : "=a" (_eax), "=S" (_ebx), "=c" (_ecx), "=d" (_edx) \
1264
                         : "0" (level));
1265
#elif defined(__i386__)
1266
        __asm__ ("push %%ebx\n\tcpuid\n\tpop %%ebx\n"                 \
1267
                         : "=a" (_eax), "=S" (_ebx), "=c" (_ecx), "=d" (_edx) \
1268
                         : "0" (level));
1269
#else
1270
        _eax=_ebx=_ecx=_edx=0;
1271
#endif
1272
    
1273
        /* Features:                                                                                                       
1274
         *                                                                                                                 
1275
         * SSE      Bit 25 of edx should be set                                                                            
1276
         * SSE2     Bit 26 of edx should be set                                                                            
1277
         * SSE3     Bit  0 of ecx should be set                                                                            
1278
         * SSE4.1   Bit 19 of ecx should be set                                                                            
1279
         */
1280
        status =  (_edx & (1 << 26)) != 0;
1281
    
1282
#else
1283
        int status = 0;
1284
#endif
1285
        /* Return SSE2 status */
1286
        return status;
1287
}
1288

    
1289

    
1290

    
1291

    
1292
void init_forcerec(FILE *fp,
1293
                   const output_env_t oenv,
1294
                   t_forcerec *fr,
1295
                   t_fcdata   *fcd,
1296
                   const t_inputrec *ir,
1297
                   const gmx_mtop_t *mtop,
1298
                   const t_commrec  *cr,
1299
                   matrix     box,
1300
                   gmx_bool       bMolEpot,
1301
                   const char *tabfn,
1302
                   const char *tabpfn,
1303
                   const char *tabbfn,
1304
                   gmx_bool       bNoSolvOpt,
1305
                   real       print_force)
1306
{
1307
    int     i,j,m,natoms,ngrp,negp_pp,negptable,egi,egj;
1308
    real    rtab;
1309
    char    *env;
1310
    double  dbl;
1311
    rvec    box_size;
1312
    const t_block *cgs;
1313
    gmx_bool    bGenericKernelOnly;
1314
    gmx_bool    bTab,bSep14tab,bNormalnblists;
1315
    t_nblists *nbl;
1316
    int     *nm_ind,egp_flags;
1317
    
1318
    fr->bDomDec = DOMAINDECOMP(cr);
1319

    
1320
    natoms = mtop->natoms;
1321

    
1322
    if (check_box(ir->ePBC,box))
1323
    {
1324
        gmx_fatal(FARGS,check_box(ir->ePBC,box));
1325
    }
1326
    
1327
    /* Test particle insertion ? */
1328
    if (EI_TPI(ir->eI)) {
1329
        /* Set to the size of the molecule to be inserted (the last one) */
1330
        /* Because of old style topologies, we have to use the last cg
1331
         * instead of the last molecule type.
1332
         */
1333
        cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
1334
        fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
1335
        if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1]) {
1336
            gmx_fatal(FARGS,"The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1337
        }
1338
    } else {
1339
        fr->n_tpi = 0;
1340
    }
1341
    
1342
    /* Copy the user determined parameters */
1343
    fr->userint1 = ir->userint1;
1344
    fr->userint2 = ir->userint2;
1345
    fr->userint3 = ir->userint3;
1346
    fr->userint4 = ir->userint4;
1347
    fr->userreal1 = ir->userreal1;
1348
    fr->userreal2 = ir->userreal2;
1349
    fr->userreal3 = ir->userreal3;
1350
    fr->userreal4 = ir->userreal4;
1351
    
1352
    /* Shell stuff */
1353
    fr->fc_stepsize = ir->fc_stepsize;
1354
    
1355
    /* Free energy */
1356
    fr->efep          = ir->efep;
1357
    fr->sc_alpha      = ir->sc_alpha;
1358
    fr->sc_power      = ir->sc_power;
1359
    fr->sc_sigma6_def = pow(ir->sc_sigma,6);
1360
    fr->sc_sigma6_min = pow(ir->sc_sigma_min,6);
1361
    env = getenv("GMX_SCSIGMA_MIN");
1362
    if (env != NULL)
1363
    {
1364
        dbl = 0;
1365
        sscanf(env,"%lf",&dbl);
1366
        fr->sc_sigma6_min = pow(dbl,6);
1367
        if (fp)
1368
        {
1369
            fprintf(fp,"Setting the minimum soft core sigma to %g nm\n",dbl);
1370
        }
1371
    }
1372

    
1373
    bGenericKernelOnly = FALSE;
1374
    if (getenv("GMX_NB_GENERIC") != NULL)
1375
    {
1376
        if (fp != NULL)
1377
        {
1378
            fprintf(fp,
1379
                    "Found environment variable GMX_NB_GENERIC.\n"
1380
                    "Disabling interaction-specific nonbonded kernels.\n\n");
1381
        }
1382
        bGenericKernelOnly = TRUE;
1383
        bNoSolvOpt         = TRUE;
1384
    }
1385
    
1386
    fr->UseOptimizedKernels = (getenv("GMX_NOOPTIMIZEDKERNELS") == NULL);
1387
    if(fp && fr->UseOptimizedKernels==FALSE)
1388
    {
1389
        fprintf(fp,
1390
                "\nFound environment variable GMX_NOOPTIMIZEDKERNELS.\n"
1391
                "Disabling SSE/SSE2/Altivec/ia64/Power6/Bluegene specific kernels.\n\n");
1392
    }    
1393

    
1394
#if ( defined(GMX_IA32_SSE2) || defined(GMX_X86_64_SSE2) || defined(GMX_IA32_SSE) || defined(GMX_X86_64_SSE)|| defined(GMX_SSE2) )
1395
    if( forcerec_check_sse2() == 0 )
1396
    {
1397
        fr->UseOptimizedKernels = FALSE;
1398
    }
1399
#endif
1400
    
1401
    /* Check if we can/should do all-vs-all kernels */
1402
    fr->bAllvsAll       = can_use_allvsall(ir,mtop,FALSE,NULL,NULL);
1403
    fr->AllvsAll_work   = NULL;
1404
    fr->AllvsAll_workgb = NULL;
1405

    
1406
    
1407
    
1408
    /* Neighbour searching stuff */
1409
    fr->bGrid      = (ir->ns_type == ensGRID);
1410
    fr->ePBC       = ir->ePBC;
1411
    fr->bMolPBC    = ir->bPeriodicMols;
1412
    fr->rc_scaling = ir->refcoord_scaling;
1413
    copy_rvec(ir->posres_com,fr->posres_com);
1414
    copy_rvec(ir->posres_comB,fr->posres_comB);
1415
    fr->rlist      = cutoff_inf(ir->rlist);
1416
    fr->rlistlong  = cutoff_inf(ir->rlistlong);
1417
    fr->eeltype    = ir->coulombtype;
1418
    fr->vdwtype    = ir->vdwtype;
1419
    
1420
    fr->bTwinRange = fr->rlistlong > fr->rlist;
1421
    fr->bEwald     = (EEL_PME(fr->eeltype) || fr->eeltype==eelEWALD);
1422
    
1423
    fr->reppow     = mtop->ffparams.reppow;
1424
    fr->bvdwtab    = (fr->vdwtype != evdwCUT ||
1425
                      !gmx_within_tol(fr->reppow,12.0,10*GMX_DOUBLE_EPS));
1426
    fr->bcoultab   = (!(fr->eeltype == eelCUT || EEL_RF(fr->eeltype)) ||
1427
                      fr->eeltype == eelRF_ZERO);
1428
    
1429
    if (getenv("GMX_FORCE_TABLES"))
1430
    {
1431
        fr->bvdwtab  = TRUE;
1432
        fr->bcoultab = TRUE;
1433
    }
1434
    
1435
    if (fp) {
1436
        fprintf(fp,"Table routines are used for coulomb: %s\n",bool_names[fr->bcoultab]);
1437
        fprintf(fp,"Table routines are used for vdw:     %s\n",bool_names[fr->bvdwtab ]);
1438
    }
1439
    
1440
    /* Tables are used for direct ewald sum */
1441
    if(fr->bEwald)
1442
    {
1443
        if (EEL_PME(ir->coulombtype))
1444
        {
1445
            if (fp)
1446
                fprintf(fp,"Will do PME sum in reciprocal space.\n");
1447
            please_cite(fp,"Essman95a");
1448
            
1449
            if (ir->ewald_geometry == eewg3DC)
1450
            {
1451
                if (fp)
1452
                {
1453
                    fprintf(fp,"Using the Ewald3DC correction for systems with a slab geometry.\n");
1454
                }
1455
                please_cite(fp,"In-Chul99a");
1456
            }
1457
        }
1458
        fr->ewaldcoeff=calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
1459
        init_ewald_tab(&(fr->ewald_table), cr, ir, fp);
1460
        if (fp)
1461
        {
1462
            fprintf(fp,"Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1463
                    1/fr->ewaldcoeff);
1464
        }
1465
    }
1466
    
1467
    /* Electrostatics */
1468
    fr->epsilon_r  = ir->epsilon_r;
1469
    fr->epsilon_rf = ir->epsilon_rf;
1470
    fr->fudgeQQ    = mtop->ffparams.fudgeQQ;
1471
    fr->rcoulomb_switch = ir->rcoulomb_switch;
1472
    fr->rcoulomb        = cutoff_inf(ir->rcoulomb);
1473
    
1474
    /* Parameters for generalized RF */
1475
    fr->zsquare = 0.0;
1476
    fr->temp    = 0.0;
1477
    
1478
    if (fr->eeltype == eelGRF)
1479
    {
1480
        init_generalized_rf(fp,mtop,ir,fr);
1481
    }
1482
    else if (EEL_FULL(fr->eeltype) || (fr->eeltype == eelSHIFT) || 
1483
             (fr->eeltype == eelUSER) || (fr->eeltype == eelSWITCH))
1484
    {
1485
        /* We must use the long range cut-off for neighboursearching...
1486
         * An extra range of e.g. 0.1 nm (half the size of a charge group)
1487
         * is necessary for neighboursearching. This allows diffusion 
1488
         * into the cut-off range (between neighborlist updates), 
1489
         * and gives more accurate forces because all atoms within the short-range
1490
         * cut-off rc must be taken into account, while the ns criterium takes
1491
         * only those with the center of geometry within the cut-off.
1492
         * (therefore we have to add half the size of a charge group, plus
1493
         * something to account for diffusion if we have nstlist > 1)
1494
         */
1495
        for(m=0; (m<DIM); m++)
1496
            box_size[m]=box[m][m];
1497
        
1498
        if (fr->eeltype == eelPPPM && fr->phi == NULL)
1499
            snew(fr->phi,natoms);
1500
        
1501
        if ((fr->eeltype==eelPPPM) || (fr->eeltype==eelPOISSON) || 
1502
            (fr->eeltype == eelSHIFT && fr->rcoulomb > fr->rcoulomb_switch))
1503
            set_shift_consts(fp,fr->rcoulomb_switch,fr->rcoulomb,box_size,fr);
1504
    }
1505
    
1506
    fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) ||
1507
                       gmx_mtop_ftype_count(mtop,F_POSRES) > 0 ||
1508
                       IR_ELEC_FIELD(*ir));
1509
    
1510
    /* Mask that says whether or not this NBF list should be computed */
1511
    /*  if (fr->bMask == NULL) {
1512
        ngrp = ir->opts.ngener*ir->opts.ngener;
1513
        snew(fr->bMask,ngrp);*/
1514
    /* Defaults to always */
1515
    /*    for(i=0; (i<ngrp); i++)
1516
          fr->bMask[i] = TRUE;
1517
          }*/
1518
    
1519
    if (ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr)) {
1520
        /* Count the total number of charge groups */
1521
        fr->cg_nalloc = ncg_mtop(mtop);
1522
        srenew(fr->cg_cm,fr->cg_nalloc);
1523
    }
1524
    if (fr->shift_vec == NULL)
1525
        snew(fr->shift_vec,SHIFTS);
1526
    
1527
    if (fr->fshift == NULL)
1528
        snew(fr->fshift,SHIFTS);
1529
    
1530
    if (fr->nbfp == NULL) {
1531
        fr->ntype = mtop->ffparams.atnr;
1532
        fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1533
        fr->nbfp  = mk_nbfp(&mtop->ffparams,fr->bBHAM);
1534
    }
1535
    
1536
    /* Copy the energy group exclusions */
1537
    fr->egp_flags = ir->opts.egp_flags;
1538
    
1539
    /* Van der Waals stuff */
1540
    fr->rvdw        = cutoff_inf(ir->rvdw);
1541
    fr->rvdw_switch = ir->rvdw_switch;
1542
    if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM) {
1543
        if (fr->rvdw_switch >= fr->rvdw)
1544
            gmx_fatal(FARGS,"rvdw_switch (%f) must be < rvdw (%f)",
1545
                      fr->rvdw_switch,fr->rvdw);
1546
        if (fp)
1547
            fprintf(fp,"Using %s Lennard-Jones, switch between %g and %g nm\n",
1548
                    (fr->eeltype==eelSWITCH) ? "switched":"shifted",
1549
                    fr->rvdw_switch,fr->rvdw);
1550
    } 
1551
    
1552
    if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
1553
        gmx_fatal(FARGS,"Switch/shift interaction not supported with Buckingham");
1554
    
1555
    if (fp)
1556
        fprintf(fp,"Cut-off's:   NS: %g   Coulomb: %g   %s: %g\n",
1557
                fr->rlist,fr->rcoulomb,fr->bBHAM ? "BHAM":"LJ",fr->rvdw);
1558
    
1559
    fr->eDispCorr = ir->eDispCorr;
1560
    if (ir->eDispCorr != edispcNO)
1561
    {
1562
        set_avcsixtwelve(fp,fr,mtop);
1563
    }
1564
    
1565
    if (fr->bBHAM)
1566
    {
1567
        set_bham_b_max(fp,fr,mtop);
1568
    }
1569

    
1570
    fr->bGB = (ir->implicit_solvent == eisGBSA);
1571
        fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
1572

    
1573
    /* Copy the GBSA data (radius, volume and surftens for each
1574
     * atomtype) from the topology atomtype section to forcerec.
1575
     */
1576
    snew(fr->atype_radius,fr->ntype);
1577
    snew(fr->atype_vol,fr->ntype);
1578
    snew(fr->atype_surftens,fr->ntype);
1579
    snew(fr->atype_gb_radius,fr->ntype);
1580
    snew(fr->atype_S_hct,fr->ntype);
1581

    
1582
    if (mtop->atomtypes.nr > 0)
1583
    {
1584
        for(i=0;i<fr->ntype;i++)
1585
            fr->atype_radius[i] =mtop->atomtypes.radius[i];
1586
        for(i=0;i<fr->ntype;i++)
1587
            fr->atype_vol[i] = mtop->atomtypes.vol[i];
1588
        for(i=0;i<fr->ntype;i++)
1589
            fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
1590
        for(i=0;i<fr->ntype;i++)
1591
            fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
1592
        for(i=0;i<fr->ntype;i++)
1593
            fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
1594
    }  
1595
        
1596
        /* Generate the GB table if needed */
1597
        if(fr->bGB)
1598
        {
1599
#ifdef GMX_DOUBLE
1600
                fr->gbtabscale=2000;
1601
#else
1602
                fr->gbtabscale=500;
1603
#endif
1604
                
1605
                fr->gbtabr=100;
1606
                fr->gbtab=make_gb_table(fp,oenv,fr,tabpfn,fr->gbtabscale);
1607

    
1608
        init_gb(&fr->born,cr,fr,ir,mtop,ir->rgbradii,ir->gb_algorithm);
1609

    
1610
        /* Copy local gb data (for dd, this is done in dd_partition_system) */
1611
        if (!DOMAINDECOMP(cr))
1612
        {
1613
            make_local_gb(cr,fr->born,ir->gb_algorithm);
1614
        }
1615
    }
1616

    
1617
    /* Set the charge scaling */
1618
    if (fr->epsilon_r != 0)
1619
        fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
1620
    else
1621
        /* eps = 0 is infinite dieletric: no coulomb interactions */
1622
        fr->epsfac = 0;
1623
    
1624
    /* Reaction field constants */
1625
    if (EEL_RF(fr->eeltype))
1626
        calc_rffac(fp,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
1627
                   fr->rcoulomb,fr->temp,fr->zsquare,box,
1628
                   &fr->kappa,&fr->k_rf,&fr->c_rf);
1629
    
1630
    set_chargesum(fp,fr,mtop);
1631
    
1632
    /* if we are using LR electrostatics, and they are tabulated,
1633
     * the tables will contain modified coulomb interactions.
1634
     * Since we want to use the non-shifted ones for 1-4
1635
     * coulombic interactions, we must have an extra set of tables.
1636
     */
1637
    
1638
    /* Construct tables.
1639
     * A little unnecessary to make both vdw and coul tables sometimes,
1640
     * but what the heck... */
1641
    
1642
    bTab = fr->bcoultab || fr->bvdwtab;
1643

    
1644
    bSep14tab = ((!bTab || fr->eeltype!=eelCUT || fr->vdwtype!=evdwCUT ||
1645
                  fr->bBHAM) &&
1646
                 (gmx_mtop_ftype_count(mtop,F_LJ14) > 0 ||
1647
                  gmx_mtop_ftype_count(mtop,F_LJC14_Q) > 0 ||
1648
                  gmx_mtop_ftype_count(mtop,F_LJC_PAIRS_NB) > 0));
1649

    
1650
    negp_pp = ir->opts.ngener - ir->nwall;
1651
    negptable = 0;
1652
    if (!bTab) {
1653
        bNormalnblists = TRUE;
1654
        fr->nnblists = 1;
1655
    } else {
1656
        bNormalnblists = (ir->eDispCorr != edispcNO);
1657
        for(egi=0; egi<negp_pp; egi++) {
1658
            for(egj=egi;  egj<negp_pp; egj++) {
1659
                egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
1660
                if (!(egp_flags & EGP_EXCL)) {
1661
                    if (egp_flags & EGP_TABLE) {
1662
                        negptable++;
1663
                    } else {
1664
                        bNormalnblists = TRUE;
1665
                    }
1666
                }
1667
            }
1668
        }
1669
        if (bNormalnblists) {
1670
            fr->nnblists = negptable + 1;
1671
        } else {
1672
            fr->nnblists = negptable;
1673
        }
1674
        if (fr->nnblists > 1)
1675
            snew(fr->gid2nblists,ir->opts.ngener*ir->opts.ngener);
1676
    }
1677
    snew(fr->nblists,fr->nnblists);
1678
    
1679
    /* This code automatically gives table length tabext without cut-off's,
1680
     * in that case grompp should already have checked that we do not need
1681
     * normal tables and we only generate tables for 1-4 interactions.
1682
     */
1683
    rtab = ir->rlistlong + ir->tabext;
1684

    
1685
    if (bTab) {
1686
        /* make tables for ordinary interactions */
1687
        if (bNormalnblists) {
1688
            make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,NULL,NULL,&fr->nblists[0]);
1689
            if (!bSep14tab)
1690
                fr->tab14 = fr->nblists[0].tab;
1691
            m = 1;
1692
        } else {
1693
            m = 0;
1694
        }
1695
        if (negptable > 0) {
1696
            /* Read the special tables for certain energy group pairs */
1697
            nm_ind = mtop->groups.grps[egcENER].nm_ind;
1698
            for(egi=0; egi<negp_pp; egi++) {
1699
                for(egj=egi;  egj<negp_pp; egj++) {
1700
                    egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
1701
                    if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL)) {
1702
                        nbl = &(fr->nblists[m]);
1703
                        if (fr->nnblists > 1) {
1704
                            fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = m;
1705
                        }
1706
                        /* Read the table file with the two energy groups names appended */
1707
                        make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,
1708
                                        *mtop->groups.grpname[nm_ind[egi]],
1709
                                        *mtop->groups.grpname[nm_ind[egj]],
1710
                                        &fr->nblists[m]);
1711
                        m++;
1712
                    } else if (fr->nnblists > 1) {
1713
                        fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = 0;
1714
                    }
1715
                }
1716
            }
1717
        }
1718
    }
1719
    if (bSep14tab)
1720
    {
1721
        /* generate extra tables with plain Coulomb for 1-4 interactions only */
1722
        fr->tab14 = make_tables(fp,oenv,fr,MASTER(cr),tabpfn,rtab,
1723
                                GMX_MAKETABLES_14ONLY);
1724
    }
1725
    
1726
    /* Wall stuff */
1727
    fr->nwall = ir->nwall;
1728
    if (ir->nwall && ir->wall_type==ewtTABLE)
1729
    {
1730
        make_wall_tables(fp,oenv,ir,tabfn,&mtop->groups,fr);
1731
    }
1732
    
1733
    if (fcd && tabbfn) {
1734
        fcd->bondtab  = make_bonded_tables(fp,
1735
                                           F_TABBONDS,F_TABBONDSNC,
1736
                                           mtop,tabbfn,"b");
1737
        fcd->angletab = make_bonded_tables(fp,
1738
                                           F_TABANGLES,-1,
1739
                                           mtop,tabbfn,"a");
1740
        fcd->dihtab   = make_bonded_tables(fp,
1741
                                           F_TABDIHS,-1,
1742
                                           mtop,tabbfn,"d");
1743
    } else {
1744
        if (debug)
1745
            fprintf(debug,"No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
1746
    }
1747
    
1748
    if (ir->eDispCorr != edispcNO)
1749
    {
1750
        calc_enervirdiff(fp,ir->eDispCorr,fr);
1751
    }
1752

    
1753
    /* QM/MM initialization if requested
1754
     */
1755
    if (ir->bQMMM)
1756
    {
1757
        fprintf(stderr,"QM/MM calculation requested.\n");
1758
    }
1759
    
1760
    fr->bQMMM      = ir->bQMMM;   
1761
    fr->qr         = mk_QMMMrec();
1762
    
1763
    /* Set all the static charge group info */
1764
    fr->cginfo_mb = init_cginfo_mb(fp,mtop,fr,bNoSolvOpt,
1765
                                   &fr->bExcl_IntraCGAll_InterCGNone);
1766
    if (DOMAINDECOMP(cr)) {
1767
        fr->cginfo = NULL;
1768
    } else {
1769
        fr->cginfo = cginfo_expand(mtop->nmolblock,fr->cginfo_mb);
1770
    }
1771
    
1772
    if (!DOMAINDECOMP(cr))
1773
    {
1774
        /* When using particle decomposition, the effect of the second argument,
1775
         * which sets fr->hcg, is corrected later in do_md and init_em.
1776
         */
1777
        forcerec_set_ranges(fr,ncg_mtop(mtop),ncg_mtop(mtop),
1778
                            mtop->natoms,mtop->natoms,mtop->natoms);
1779
    }
1780
    
1781
    fr->print_force = print_force;
1782

    
1783

    
1784
    /* coarse load balancing vars */
1785
    fr->t_fnbf=0.;
1786
    fr->t_wait=0.;
1787
    fr->timesteps=0;
1788
    
1789
    /* Initialize neighbor search */
1790
    init_ns(fp,cr,&fr->ns,fr,mtop,box);
1791
    
1792
    if (cr->duty & DUTY_PP)
1793
        gmx_setup_kernels(fp,bGenericKernelOnly);
1794
}
1795

    
1796
#define pr_real(fp,r) fprintf(fp,"%s: %e\n",#r,r)
1797
#define pr_int(fp,i)  fprintf((fp),"%s: %d\n",#i,i)
1798
#define pr_bool(fp,b) fprintf((fp),"%s: %s\n",#b,bool_names[b])
1799

    
1800
void pr_forcerec(FILE *fp,t_forcerec *fr,t_commrec *cr)
1801
{
1802
  int i;
1803

    
1804
  pr_real(fp,fr->rlist);
1805
  pr_real(fp,fr->rcoulomb);
1806
  pr_real(fp,fr->fudgeQQ);
1807
  pr_bool(fp,fr->bGrid);
1808
  pr_bool(fp,fr->bTwinRange);
1809
  /*pr_int(fp,fr->cg0);
1810
    pr_int(fp,fr->hcg);*/
1811
  for(i=0; i<fr->nnblists; i++)
1812
    pr_int(fp,fr->nblists[i].tab.n);
1813
  pr_real(fp,fr->rcoulomb_switch);
1814
  pr_real(fp,fr->rcoulomb);
1815
  
1816
  fflush(fp);
1817
}