1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
49 #include "nonbonded.h"
63 #include "mtop_util.h"
65 t_forcerec
*mk_forcerec(void)
75 static void pr_nbfp(FILE *fp
,real
*nbfp
,bool bBHAM
,int atnr
)
79 for(i
=0; (i
<atnr
); i
++) {
80 for(j
=0; (j
<atnr
); j
++) {
81 fprintf(fp
,"%2d - %2d",i
,j
);
83 fprintf(fp
," a=%10g, b=%10g, c=%10g\n",BHAMA(nbfp
,atnr
,i
,j
),
84 BHAMB(nbfp
,atnr
,i
,j
),BHAMC(nbfp
,atnr
,i
,j
));
86 fprintf(fp
," c6=%10g, c12=%10g\n",C6(nbfp
,atnr
,i
,j
),
93 static real
*mk_nbfp(const gmx_ffparams_t
*idef
,bool bBHAM
)
100 snew(nbfp
,3*atnr
*atnr
);
101 for(i
=k
=0; (i
<atnr
); i
++) {
102 for(j
=0; (j
<atnr
); j
++,k
++) {
103 BHAMA(nbfp
,atnr
,i
,j
) = idef
->iparams
[k
].bham
.a
;
104 BHAMB(nbfp
,atnr
,i
,j
) = idef
->iparams
[k
].bham
.b
;
105 BHAMC(nbfp
,atnr
,i
,j
) = idef
->iparams
[k
].bham
.c
;
110 snew(nbfp
,2*atnr
*atnr
);
111 for(i
=k
=0; (i
<atnr
); i
++) {
112 for(j
=0; (j
<atnr
); j
++,k
++) {
113 C6(nbfp
,atnr
,i
,j
) = idef
->iparams
[k
].lj
.c6
;
114 C12(nbfp
,atnr
,i
,j
) = idef
->iparams
[k
].lj
.c12
;
121 /* This routine sets fr->solvent_opt to the most common solvent in the
122 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
123 * the fr->solvent_type array with the correct type (or esolNO).
125 * Charge groups that fulfill the conditions but are not identical to the
126 * most common one will be marked as esolNO in the solvent_type array.
128 * TIP3p is identical to SPC for these purposes, so we call it
129 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
131 * NOTE: QM particle should not
132 * become an optimized solvent. Not even if there is only one charge
142 } solvent_parameters_t
;
145 check_solvent_cg(const gmx_moltype_t
*molt
,
148 const unsigned char *qm_grpnr
,
149 const t_grps
*qm_grps
,
151 int *n_solvent_parameters
,
152 solvent_parameters_t
**solvent_parameters_p
,
156 const t_blocka
* excl
;
167 solvent_parameters_t
*solvent_parameters
;
169 /* We use a list with parameters for each solvent type.
170 * Every time we discover a new molecule that fulfills the basic
171 * conditions for a solvent we compare with the previous entries
172 * in these lists. If the parameters are the same we just increment
173 * the counter for that type, and otherwise we create a new type
174 * based on the current molecule.
176 * Once we've finished going through all molecules we check which
177 * solvent is most common, and mark all those molecules while we
178 * clear the flag on all others.
181 solvent_parameters
= *solvent_parameters_p
;
183 /* Mark the cg first as non optimized */
186 /* Check if this cg has no exclusions with atoms in other charge groups
187 * and all atoms inside the charge group excluded.
188 * We only have 3 or 4 atom solvent loops.
190 if (GET_CGINFO_EXCL_INTER(cginfo
) ||
191 !GET_CGINFO_EXCL_INTRA(cginfo
))
196 /* Get the indices of the first atom in this charge group */
197 j0
= molt
->cgs
.index
[cg0
];
198 j1
= molt
->cgs
.index
[cg0
+1];
200 /* Number of atoms in our molecule */
205 "Moltype '%s': there are %d atoms in this charge group\n",
209 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
217 /* Check if we are doing QM on this group */
219 if (qm_grpnr
!= NULL
)
221 for(j
=j0
; j
<j1
&& !qm
; j
++)
223 qm
= (qm_grpnr
[j
] < qm_grps
->nr
- 1);
226 /* Cannot use solvent optimization with QM */
232 atom
= molt
->atoms
.atom
;
234 /* Still looks like a solvent, time to check parameters */
236 /* If it is perturbed (free energy) we can't use the solvent loops,
237 * so then we just skip to the next molecule.
241 for(j
=j0
; j
<j1
&& !perturbed
; j
++)
243 perturbed
= PERTURBED(atom
[j
]);
251 /* Now it's only a question if the VdW and charge parameters
252 * are OK. Before doing the check we compare and see if they are
253 * identical to a possible previous solvent type.
254 * First we assign the current types and charges.
258 tmp_vdwtype
[j
] = atom
[j0
+j
].type
;
259 tmp_charge
[j
] = atom
[j0
+j
].q
;
262 /* Does it match any previous solvent type? */
263 for(k
=0 ; k
<*n_solvent_parameters
; k
++)
268 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
269 if( (solvent_parameters
[k
].model
==esolSPC
&& nj
!=3) ||
270 (solvent_parameters
[k
].model
==esolTIP4P
&& nj
!=4) )
273 /* Check that types & charges match for all atoms in molecule */
274 for(j
=0 ; j
<nj
&& match
==TRUE
; j
++)
276 if (tmp_vdwtype
[j
] != solvent_parameters
[k
].vdwtype
[j
])
280 if(tmp_charge
[j
] != solvent_parameters
[k
].charge
[j
])
287 /* Congratulations! We have a matched solvent.
288 * Flag it with this type for later processing.
291 solvent_parameters
[k
].count
+= nmol
;
293 /* We are done with this charge group */
298 /* If we get here, we have a tentative new solvent type.
299 * Before we add it we must check that it fulfills the requirements
300 * of the solvent optimized loops. First determine which atoms have
306 tjA
= tmp_vdwtype
[j
];
308 /* Go through all other tpes and see if any have non-zero
309 * VdW parameters when combined with this one.
311 for(k
=0; k
<fr
->ntype
&& (has_vdw
[j
]==FALSE
); k
++)
313 /* We already checked that the atoms weren't perturbed,
314 * so we only need to check state A now.
318 has_vdw
[j
] = (has_vdw
[j
] ||
319 (BHAMA(fr
->nbfp
,fr
->ntype
,tjA
,k
) != 0.0) ||
320 (BHAMB(fr
->nbfp
,fr
->ntype
,tjA
,k
) != 0.0) ||
321 (BHAMC(fr
->nbfp
,fr
->ntype
,tjA
,k
) != 0.0));
326 has_vdw
[j
] = (has_vdw
[j
] ||
327 (C6(fr
->nbfp
,fr
->ntype
,tjA
,k
) != 0.0) ||
328 (C12(fr
->nbfp
,fr
->ntype
,tjA
,k
) != 0.0));
333 /* Now we know all we need to make the final check and assignment. */
337 * For this we require thatn all atoms have charge,
338 * the charges on atom 2 & 3 should be the same, and only
339 * atom 1 should have VdW.
341 if (has_vdw
[0] == TRUE
&&
342 has_vdw
[1] == FALSE
&&
343 has_vdw
[2] == FALSE
&&
344 tmp_charge
[0] != 0 &&
345 tmp_charge
[1] != 0 &&
346 tmp_charge
[2] == tmp_charge
[1])
348 srenew(solvent_parameters
,*n_solvent_parameters
+1);
349 solvent_parameters
[*n_solvent_parameters
].model
= esolSPC
;
350 solvent_parameters
[*n_solvent_parameters
].count
= nmol
;
353 solvent_parameters
[*n_solvent_parameters
].vdwtype
[k
] = tmp_vdwtype
[k
];
354 solvent_parameters
[*n_solvent_parameters
].charge
[k
] = tmp_charge
[k
];
357 *cg_sp
= *n_solvent_parameters
;
358 (*n_solvent_parameters
)++;
363 /* Or could it be a TIP4P?
364 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
365 * Only atom 1 should have VdW.
367 if(has_vdw
[0] == TRUE
&&
368 has_vdw
[1] == FALSE
&&
369 has_vdw
[2] == FALSE
&&
370 has_vdw
[3] == FALSE
&&
371 tmp_charge
[0] == 0 &&
372 tmp_charge
[1] != 0 &&
373 tmp_charge
[2] == tmp_charge
[1] &&
376 srenew(solvent_parameters
,*n_solvent_parameters
+1);
377 solvent_parameters
[*n_solvent_parameters
].model
= esolTIP4P
;
378 solvent_parameters
[*n_solvent_parameters
].count
= nmol
;
381 solvent_parameters
[*n_solvent_parameters
].vdwtype
[k
] = tmp_vdwtype
[k
];
382 solvent_parameters
[*n_solvent_parameters
].charge
[k
] = tmp_charge
[k
];
385 *cg_sp
= *n_solvent_parameters
;
386 (*n_solvent_parameters
)++;
390 *solvent_parameters_p
= solvent_parameters
;
394 check_solvent(FILE * fp
,
395 const gmx_mtop_t
* mtop
,
397 cginfo_mb_t
*cginfo_mb
)
400 const t_block
* mols
;
401 const gmx_moltype_t
*molt
;
402 int mb
,mol
,cg_mol
,at_offset
,cg_offset
,am
,cgm
,i
,nmol_ch
,nmol
;
403 int n_solvent_parameters
;
404 solvent_parameters_t
*solvent_parameters
;
410 fprintf(debug
,"Going to determine what solvent types we have.\n");
415 n_solvent_parameters
= 0;
416 solvent_parameters
= NULL
;
417 /* Allocate temporary array for solvent type */
418 snew(cg_sp
,mtop
->nmolblock
);
422 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
424 molt
= &mtop
->moltype
[mtop
->molblock
[mb
].type
];
426 /* Here we have to loop over all individual molecules
427 * because we need to check for QMMM particles.
429 snew(cg_sp
[mb
],cginfo_mb
[mb
].cg_mod
);
430 nmol_ch
= cginfo_mb
[mb
].cg_mod
/cgs
->nr
;
431 nmol
= mtop
->molblock
[mb
].nmol
/nmol_ch
;
432 for(mol
=0; mol
<nmol_ch
; mol
++)
435 am
= mol
*cgs
->index
[cgs
->nr
];
436 for(cg_mol
=0; cg_mol
<cgs
->nr
; cg_mol
++)
438 check_solvent_cg(molt
,cg_mol
,nmol
,
439 mtop
->groups
.grpnr
[egcQMMM
] ?
440 mtop
->groups
.grpnr
[egcQMMM
]+at_offset
+am
: 0,
441 &mtop
->groups
.grps
[egcQMMM
],
443 &n_solvent_parameters
,&solvent_parameters
,
444 cginfo_mb
[mb
].cginfo
[cgm
+cg_mol
],
445 &cg_sp
[mb
][cgm
+cg_mol
]);
448 cg_offset
+= cgs
->nr
;
449 at_offset
+= cgs
->index
[cgs
->nr
];
452 /* Puh! We finished going through all charge groups.
453 * Now find the most common solvent model.
456 /* Most common solvent this far */
458 for(i
=0;i
<n_solvent_parameters
;i
++)
461 solvent_parameters
[i
].count
> solvent_parameters
[bestsp
].count
)
469 bestsol
= solvent_parameters
[bestsp
].model
;
476 #ifdef DISABLE_WATER_NLIST
481 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
483 cgs
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].cgs
;
484 nmol
= (mtop
->molblock
[mb
].nmol
*cgs
->nr
)/cginfo_mb
[mb
].cg_mod
;
485 for(i
=0; i
<cginfo_mb
[mb
].cg_mod
; i
++)
487 if (cg_sp
[mb
][i
] == bestsp
)
489 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[i
],bestsol
);
494 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[i
],esolNO
);
501 if (bestsol
!= esolNO
&& fp
!=NULL
)
503 fprintf(fp
,"\nEnabling %s-like water optimization for %d molecules.\n\n",
505 solvent_parameters
[bestsp
].count
);
508 sfree(solvent_parameters
);
509 fr
->solvent_opt
= bestsol
;
512 static cginfo_mb_t
*init_cginfo_mb(FILE *fplog
,const gmx_mtop_t
*mtop
,
513 t_forcerec
*fr
,bool bNoSolvOpt
)
516 const t_blocka
*excl
;
517 const gmx_moltype_t
*molt
;
518 const gmx_molblock_t
*molb
;
519 cginfo_mb_t
*cginfo_mb
;
521 int cg_offset
,a_offset
,cgm
,am
;
522 int mb
,m
,ncg_tot
,cg
,a0
,a1
,gid
,ai
,j
,aj
,excl_nalloc
;
523 bool bId
,*bExcl
,bExclIntraAll
,bExclInter
;
525 ncg_tot
= ncg_mtop(mtop
);
526 snew(cginfo_mb
,mtop
->nmolblock
);
529 snew(bExcl
,excl_nalloc
);
532 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
534 molb
= &mtop
->molblock
[mb
];
535 molt
= &mtop
->moltype
[molb
->type
];
539 /* Check if the cginfo is identical for all molecules in this block.
540 * If so, we only need an array of the size of one molecule.
541 * Otherwise we make an array of #mol times #cgs per molecule.
545 for(m
=0; m
<molb
->nmol
; m
++)
547 am
= m
*cgs
->index
[cgs
->nr
];
548 for(cg
=0; cg
<cgs
->nr
; cg
++)
551 a1
= cgs
->index
[cg
+1];
552 if (ggrpnr(&mtop
->groups
,egcENER
,a_offset
+am
+a0
) !=
553 ggrpnr(&mtop
->groups
,egcENER
,a_offset
+a0
))
557 if (mtop
->groups
.grpnr
[egcQMMM
] != NULL
)
559 for(ai
=a0
; ai
<a1
; ai
++)
561 if (mtop
->groups
.grpnr
[egcQMMM
][a_offset
+am
+ai
] !=
562 mtop
->groups
.grpnr
[egcQMMM
][a_offset
+ai
])
571 cginfo_mb
[mb
].cg_start
= cg_offset
;
572 cginfo_mb
[mb
].cg_end
= cg_offset
+ molb
->nmol
*cgs
->nr
;
573 cginfo_mb
[mb
].cg_mod
= (bId
? 1 : molb
->nmol
)*cgs
->nr
;
574 snew(cginfo_mb
[mb
].cginfo
,cginfo_mb
[mb
].cg_mod
);
575 cginfo
= cginfo_mb
[mb
].cginfo
;
577 for(m
=0; m
<(bId
? 1 : molb
->nmol
); m
++)
580 am
= m
*cgs
->index
[cgs
->nr
];
581 for(cg
=0; cg
<cgs
->nr
; cg
++)
584 a1
= cgs
->index
[cg
+1];
586 /* Store the energy group in cginfo */
587 gid
= ggrpnr(&mtop
->groups
,egcENER
,a_offset
+am
+a0
);
588 SET_CGINFO_GID(cginfo
[cgm
+cg
],gid
);
590 /* Check the intra/inter charge group exclusions */
591 if (a1
-a0
> excl_nalloc
) {
592 excl_nalloc
= a1
- a0
;
593 srenew(bExcl
,excl_nalloc
);
595 /* bExclIntraAll: all intra cg interactions excluded
596 * bExclInter: any inter cg interactions excluded
598 bExclIntraAll
= TRUE
;
600 for(ai
=a0
; ai
<a1
; ai
++) {
601 /* Clear the exclusion list for atom ai */
602 for(aj
=a0
; aj
<a1
; aj
++) {
603 bExcl
[aj
-a0
] = FALSE
;
605 /* Loop over all the exclusions of atom ai */
606 for(j
=excl
->index
[ai
]; j
<excl
->index
[ai
+1]; j
++)
609 if (aj
< a0
|| aj
>= a1
)
618 /* Check if ai excludes a0 to a1 */
619 for(aj
=a0
; aj
<a1
; aj
++)
623 bExclIntraAll
= FALSE
;
629 SET_CGINFO_EXCL_INTRA(cginfo
[cgm
+cg
]);
633 SET_CGINFO_EXCL_INTER(cginfo
[cgm
+cg
]);
635 if (a1
- a0
> MAX_CHARGEGROUP_SIZE
)
637 /* The size in cginfo is currently only read with DD */
638 gmx_fatal(FARGS
,"A charge group has size %d which is larger than the limit of %d atoms",a1
-a0
,MAX_CHARGEGROUP_SIZE
);
640 SET_CGINFO_NATOMS(cginfo
[cgm
+cg
],a1
-a0
);
643 cg_offset
+= molb
->nmol
*cgs
->nr
;
644 a_offset
+= molb
->nmol
*cgs
->index
[cgs
->nr
];
648 /* the solvent optimizer is called after the QM is initialized,
649 * because we don't want to have the QM subsystemto become an
653 check_solvent(fplog
,mtop
,fr
,cginfo_mb
);
655 if (getenv("GMX_NO_SOLV_OPT"))
659 fprintf(fplog
,"Found environment variable GMX_NO_SOLV_OPT.\n"
660 "Disabling all solvent optimization\n");
662 fr
->solvent_opt
= esolNO
;
666 fr
->solvent_opt
= esolNO
;
668 if (!fr
->solvent_opt
)
670 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
672 for(cg
=0; cg
<cginfo_mb
[mb
].cg_mod
; cg
++)
674 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[cg
],esolNO
);
682 static int *cginfo_expand(int nmb
,cginfo_mb_t
*cgi_mb
)
687 ncg
= cgi_mb
[nmb
-1].cg_end
;
690 for(cg
=0; cg
<ncg
; cg
++)
692 while (cg
>= cgi_mb
[mb
].cg_end
)
697 cgi_mb
[mb
].cginfo
[(cg
- cgi_mb
[mb
].cg_start
) % cgi_mb
[mb
].cg_mod
];
703 static void set_chargesum(FILE *log
,t_forcerec
*fr
,const gmx_mtop_t
*mtop
)
707 const t_atoms
*atoms
;
710 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
712 nmol
= mtop
->molblock
[mb
].nmol
;
713 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
714 for(i
=0; i
<atoms
->nr
; i
++)
716 qsum
+= nmol
*atoms
->atom
[i
].q
;
720 if (fr
->efep
!= efepNO
)
723 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
725 nmol
= mtop
->molblock
[mb
].nmol
;
726 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
727 for(i
=0; i
<atoms
->nr
; i
++)
729 qsum
+= nmol
*atoms
->atom
[i
].qB
;
736 fr
->qsum
[1] = fr
->qsum
[0];
739 if (fr
->efep
== efepNO
)
740 fprintf(log
,"System total charge: %.3f\n",fr
->qsum
[0]);
742 fprintf(log
,"System total charge, top. A: %.3f top. B: %.3f\n",
743 fr
->qsum
[0],fr
->qsum
[1]);
747 void update_forcerec(FILE *log
,t_forcerec
*fr
,matrix box
)
749 if (fr
->eeltype
== eelGRF
)
751 calc_rffac(NULL
,fr
->eeltype
,fr
->epsilon_r
,fr
->epsilon_rf
,
752 fr
->rcoulomb
,fr
->temp
,fr
->zsquare
,box
,
753 &fr
->kappa
,&fr
->k_rf
,&fr
->c_rf
);
757 void set_avcsixtwelve(FILE *fplog
,t_forcerec
*fr
,const gmx_mtop_t
*mtop
)
759 const t_atoms
*atoms
;
760 const t_blocka
*excl
;
761 int mb
,nmol
,nmolc
,i
,j
,tpi
,tpj
,j1
,j2
,k
,n
,nexcl
,q
;
762 #if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)
763 long long int npair
,npair_ij
,tmpi
,tmpj
;
765 double npair
, npair_ij
,tmpi
,tmpj
;
776 for(q
=0; q
<(fr
->efep
==efepNO
? 1 : 2); q
++) {
782 /* Count the types so we avoid natoms^2 operations */
784 for(mb
=0; mb
<mtop
->nmolblock
; mb
++) {
785 nmol
= mtop
->molblock
[mb
].nmol
;
786 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
787 for(i
=0; i
<atoms
->nr
; i
++) {
790 tpi
= atoms
->atom
[i
].type
;
794 tpi
= atoms
->atom
[i
].typeB
;
796 typecount
[tpi
] += nmol
;
799 for(tpi
=0; tpi
<ntp
; tpi
++) {
800 for(tpj
=tpi
; tpj
<ntp
; tpj
++) {
801 tmpi
= typecount
[tpi
];
802 tmpj
= typecount
[tpj
];
805 npair_ij
= tmpi
*tmpj
;
809 npair_ij
= tmpi
*(tmpi
- 1)/2;
812 csix
+= npair_ij
*BHAMC(nbfp
,ntp
,tpi
,tpj
);
814 csix
+= npair_ij
* C6(nbfp
,ntp
,tpi
,tpj
);
815 ctwelve
+= npair_ij
* C12(nbfp
,ntp
,tpi
,tpj
);
821 /* Subtract the excluded pairs.
822 * The main reason for substracting exclusions is that in some cases
823 * some combinations might never occur and the parameters could have
824 * any value. These unused values should not influence the dispersion
827 for(mb
=0; mb
<mtop
->nmolblock
; mb
++) {
828 nmol
= mtop
->molblock
[mb
].nmol
;
829 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
830 excl
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].excls
;
831 for(i
=0; (i
<atoms
->nr
); i
++) {
834 tpi
= atoms
->atom
[i
].type
;
838 tpi
= atoms
->atom
[i
].typeB
;
841 j2
= excl
->index
[i
+1];
842 for(j
=j1
; j
<j2
; j
++) {
848 tpj
= atoms
->atom
[k
].type
;
852 tpj
= atoms
->atom
[k
].typeB
;
855 csix
-= nmol
*BHAMC(nbfp
,ntp
,tpi
,tpj
);
857 csix
-= nmol
*C6 (nbfp
,ntp
,tpi
,tpj
);
858 ctwelve
-= nmol
*C12(nbfp
,ntp
,tpi
,tpj
);
866 /* Only correct for the interaction of the test particle
867 * with the rest of the system.
869 atoms
= &mtop
->moltype
[mtop
->molblock
[mtop
->nmolblock
-1].type
].atoms
;
872 tpi
= atoms
->atom
[atoms
->nr
-1].type
;
876 tpi
= atoms
->atom
[atoms
->nr
-1].typeB
;
879 for(mb
=0; mb
<mtop
->nmolblock
; mb
++) {
880 nmol
= mtop
->molblock
[mb
].nmol
;
881 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
882 for(j
=0; j
<atoms
->nr
; j
++) {
884 /* Remove the interaction of the test charge group
887 if (mb
== mtop
->nmolblock
-1 && j
>= atoms
->nr
- fr
->n_tpi
)
893 tpj
= atoms
->atom
[j
].type
;
897 tpj
= atoms
->atom
[j
].typeB
;
901 csix
+= nmolc
*BHAMC(nbfp
,ntp
,tpi
,tpj
);
905 csix
+= nmolc
*C6 (nbfp
,ntp
,tpi
,tpj
);
906 ctwelve
+= nmolc
*C12(nbfp
,ntp
,tpi
,tpj
);
912 if (npair
- nexcl
<= 0 && fplog
) {
913 fprintf(fplog
,"\nWARNING: There are no atom pairs for dispersion correction\n\n");
917 csix
/= npair
- nexcl
;
918 ctwelve
/= npair
- nexcl
;
921 fprintf(debug
,"Counted %d exclusions\n",nexcl
);
922 fprintf(debug
,"Average C6 parameter is: %10g\n",(double)csix
);
923 fprintf(debug
,"Average C12 parameter is: %10g\n",(double)ctwelve
);
925 fr
->avcsix
[q
] = csix
;
926 fr
->avctwelve
[q
] = ctwelve
;
930 if (fr
->eDispCorr
== edispcAllEner
||
931 fr
->eDispCorr
== edispcAllEnerPres
)
933 fprintf(fplog
,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
934 fr
->avcsix
[0],fr
->avctwelve
[0]);
938 fprintf(fplog
,"Long Range LJ corr.: <C6> %10.4e\n",fr
->avcsix
[0]);
944 static void set_bham_b_max(FILE *fplog
,t_forcerec
*fr
,
945 const gmx_mtop_t
*mtop
)
947 const t_atoms
*at1
,*at2
;
948 int mt1
,mt2
,i
,j
,tpi
,tpj
,ntypes
;
954 fprintf(fplog
,"Determining largest Buckingham b parameter for table\n");
961 for(mt1
=0; mt1
<mtop
->nmoltype
; mt1
++)
963 at1
= &mtop
->moltype
[mt1
].atoms
;
964 for(i
=0; (i
<at1
->nr
); i
++)
966 tpi
= at1
->atom
[i
].type
;
968 gmx_fatal(FARGS
,"Atomtype[%d] = %d, maximum = %d",i
,tpi
,ntypes
);
970 for(mt2
=mt1
; mt2
<mtop
->nmoltype
; mt2
++)
972 at2
= &mtop
->moltype
[mt2
].atoms
;
973 for(j
=0; (j
<at2
->nr
); j
++) {
974 tpj
= at2
->atom
[j
].type
;
977 gmx_fatal(FARGS
,"Atomtype[%d] = %d, maximum = %d",j
,tpj
,ntypes
);
979 b
= BHAMB(nbfp
,ntypes
,tpi
,tpj
);
980 if (b
> fr
->bham_b_max
)
984 if ((b
< bmin
) || (bmin
==-1))
994 fprintf(fplog
,"Buckingham b parameters, min: %g, max: %g\n",
995 bmin
,fr
->bham_b_max
);
999 static void make_nbf_tables(FILE *fp
,const output_env_t oenv
,
1000 t_forcerec
*fr
,real rtab
,
1001 const t_commrec
*cr
,
1002 const char *tabfn
,char *eg1
,char *eg2
,
1009 if (tabfn
== NULL
) {
1011 fprintf(debug
,"No table file name passed, can not read table, can not do non-bonded interactions\n");
1015 sprintf(buf
,"%s",tabfn
);
1017 /* Append the two energy group names */
1018 sprintf(buf
+ strlen(tabfn
) - strlen(ftp2ext(efXVG
)) - 1,"_%s_%s.%s",
1019 eg1
,eg2
,ftp2ext(efXVG
));
1020 nbl
->tab
= make_tables(fp
,oenv
,fr
,MASTER(cr
),buf
,rtab
,0);
1021 /* Copy the contents of the table to separate coulomb and LJ tables too,
1022 * to improve cache performance.
1025 /* For performance reasons we want
1026 * the table data to be aligned to 16-byte. This is accomplished
1027 * by allocating 16 bytes extra to a temporary pointer, and then
1028 * calculating an aligned pointer. This new pointer must not be
1029 * used in a free() call, but thankfully we're sloppy enough not
1033 /* 8 fp entries per vdw table point, n+1 points, and 16 bytes extra to align it. */
1034 p_tmp
= malloc(8*(nbl
->tab
.n
+1)*sizeof(real
)+16);
1036 /* align it - size_t has the same same as a pointer */
1037 nbl
->vdwtab
= (real
*) (((size_t) p_tmp
+ 16) & (~((size_t) 15)));
1039 /* 4 fp entries per coul table point, n+1 points, and 16 bytes extra to align it. */
1040 p_tmp
= malloc(4*(nbl
->tab
.n
+1)*sizeof(real
)+16);
1042 /* align it - size_t has the same same as a pointer */
1043 nbl
->coultab
= (real
*) (((size_t) p_tmp
+ 16) & (~((size_t) 15)));
1046 for(i
=0; i
<=nbl
->tab
.n
; i
++) {
1048 nbl
->coultab
[4*i
+j
] = nbl
->tab
.tab
[12*i
+j
];
1050 nbl
->vdwtab
[8*i
+j
] = nbl
->tab
.tab
[12*i
+4+j
];
1054 static void count_tables(int ftype1
,int ftype2
,const gmx_mtop_t
*mtop
,
1055 int *ncount
,int **count
)
1057 const gmx_moltype_t
*molt
;
1059 int mt
,ftype
,stride
,i
,j
,tabnr
;
1061 for(mt
=0; mt
<mtop
->nmoltype
; mt
++)
1063 molt
= &mtop
->moltype
[mt
];
1064 for(ftype
=0; ftype
<F_NRE
; ftype
++)
1066 if (ftype
== ftype1
|| ftype
== ftype2
) {
1067 il
= &molt
->ilist
[ftype
];
1068 stride
= 1 + NRAL(ftype
);
1069 for(i
=0; i
<il
->nr
; i
+=stride
) {
1070 tabnr
= mtop
->ffparams
.iparams
[il
->iatoms
[i
]].tab
.table
;
1072 gmx_fatal(FARGS
,"A bonded table number is smaller than 0: %d\n",tabnr
);
1073 if (tabnr
>= *ncount
) {
1074 srenew(*count
,tabnr
+1);
1075 for(j
=*ncount
; j
<tabnr
+1; j
++)
1086 static bondedtable_t
*make_bonded_tables(FILE *fplog
,
1087 int ftype1
,int ftype2
,
1088 const gmx_mtop_t
*mtop
,
1089 const char *basefn
,const char *tabext
)
1091 int i
,ncount
,*count
;
1099 count_tables(ftype1
,ftype2
,mtop
,&ncount
,&count
);
1103 for(i
=0; i
<ncount
; i
++) {
1105 sprintf(tabfn
,"%s",basefn
);
1106 sprintf(tabfn
+ strlen(basefn
) - strlen(ftp2ext(efXVG
)) - 1,"_%s%d.%s",
1107 tabext
,i
,ftp2ext(efXVG
));
1108 tab
[i
] = make_bonded_table(fplog
,tabfn
,NRAL(ftype1
)-2);
1117 void forcerec_set_ranges(t_forcerec
*fr
,
1118 int ncg_home
,int ncg_force
,
1119 int natoms_force
,int natoms_f_novirsum
)
1124 /* fr->ncg_force is unused in the standard code,
1125 * but it can be useful for modified code dealing with charge groups.
1127 fr
->ncg_force
= ncg_force
;
1128 fr
->natoms_force
= natoms_force
;
1130 if (fr
->natoms_force
> fr
->nalloc_force
)
1132 fr
->nalloc_force
= over_alloc_dd(fr
->natoms_force
);
1136 srenew(fr
->f_twin
,fr
->nalloc_force
);
1140 if (fr
->bF_NoVirSum
)
1142 fr
->f_novirsum_n
= natoms_f_novirsum
;
1143 if (fr
->f_novirsum_n
> fr
->f_novirsum_nalloc
)
1145 fr
->f_novirsum_nalloc
= over_alloc_dd(fr
->f_novirsum_n
);
1146 srenew(fr
->f_novirsum_alloc
,fr
->f_novirsum_nalloc
);
1151 fr
->f_novirsum_n
= 0;
1155 static real
cutoff_inf(real cutoff
)
1159 cutoff
= GMX_CUTOFF_INF
;
1165 void init_forcerec(FILE *fp
,
1166 const output_env_t oenv
,
1169 const t_inputrec
*ir
,
1170 const gmx_mtop_t
*mtop
,
1171 const t_commrec
*cr
,
1180 int i
,j
,m
,natoms
,ngrp
,negp_pp
,negptable
,egi
,egj
;
1184 bool bTab
,bSep14tab
,bNormalnblists
;
1186 int *nm_ind
,egp_flags
;
1188 fr
->bDomDec
= DOMAINDECOMP(cr
);
1190 natoms
= mtop
->natoms
;
1192 if (check_box(ir
->ePBC
,box
))
1194 gmx_fatal(FARGS
,check_box(ir
->ePBC
,box
));
1197 /* Test particle insertion ? */
1198 if (EI_TPI(ir
->eI
)) {
1199 /* Set to the size of the molecule to be inserted (the last one) */
1200 /* Because of old style topologies, we have to use the last cg
1201 * instead of the last molecule type.
1203 cgs
= &mtop
->moltype
[mtop
->molblock
[mtop
->nmolblock
-1].type
].cgs
;
1204 fr
->n_tpi
= cgs
->index
[cgs
->nr
] - cgs
->index
[cgs
->nr
-1];
1205 if (fr
->n_tpi
!= mtop
->mols
.index
[mtop
->mols
.nr
] - mtop
->mols
.index
[mtop
->mols
.nr
-1]) {
1206 gmx_fatal(FARGS
,"The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1212 /* Copy the user determined parameters */
1213 fr
->userint1
= ir
->userint1
;
1214 fr
->userint2
= ir
->userint2
;
1215 fr
->userint3
= ir
->userint3
;
1216 fr
->userint4
= ir
->userint4
;
1217 fr
->userreal1
= ir
->userreal1
;
1218 fr
->userreal2
= ir
->userreal2
;
1219 fr
->userreal3
= ir
->userreal3
;
1220 fr
->userreal4
= ir
->userreal4
;
1223 fr
->fc_stepsize
= ir
->fc_stepsize
;
1226 fr
->efep
= ir
->efep
;
1227 fr
->sc_alpha
= ir
->sc_alpha
;
1228 fr
->sc_power
= ir
->sc_power
;
1229 fr
->sc_sigma6
= pow(ir
->sc_sigma
,6);
1231 /* Check if we can/should do all-vs-all kernels */
1233 /* double not done yet */
1234 fr
->bAllvsAll
= FALSE
;
1236 fr
->bAllvsAll
= (ir
->rlist
==0 &&
1239 ir
->ePBC
==epbcNONE
&&
1240 ir
->vdwtype
==evdwCUT
&&
1241 ir
->coulombtype
==eelCUT
&&
1243 (ir
->implicit_solvent
== eisNO
||
1244 (ir
->implicit_solvent
==eisGBSA
&& (ir
->gb_algorithm
==egbSTILL
||
1245 ir
->gb_algorithm
==egbHCT
||
1246 ir
->gb_algorithm
==egbOBC
)))
1248 if (fr
->bAllvsAll
&& ir
->opts
.ngener
> 1)
1250 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";
1253 fprintf(stderr
,"\n%s\n",note
);
1257 fprintf(fp
,"\n%s\n",note
);
1259 fr
->bAllvsAll
= FALSE
;
1262 fr
->AllvsAll_work
= NULL
;
1263 fr
->AllvsAll_workgb
= NULL
;
1265 /* Neighbour searching stuff */
1266 fr
->bGrid
= (ir
->ns_type
== ensGRID
);
1267 fr
->ePBC
= ir
->ePBC
;
1268 fr
->bMolPBC
= ir
->bPeriodicMols
;
1269 fr
->rc_scaling
= ir
->refcoord_scaling
;
1270 copy_rvec(ir
->posres_com
,fr
->posres_com
);
1271 copy_rvec(ir
->posres_comB
,fr
->posres_comB
);
1272 fr
->rlist
= cutoff_inf(ir
->rlist
);
1273 fr
->rlistlong
= cutoff_inf(ir
->rlistlong
);
1274 fr
->eeltype
= ir
->coulombtype
;
1275 fr
->vdwtype
= ir
->vdwtype
;
1277 fr
->bTwinRange
= fr
->rlistlong
> fr
->rlist
;
1278 fr
->bEwald
= (EEL_PME(fr
->eeltype
) || fr
->eeltype
==eelEWALD
);
1280 fr
->reppow
= mtop
->ffparams
.reppow
;
1281 fr
->bvdwtab
= (fr
->vdwtype
!= evdwCUT
||
1282 !gmx_within_tol(fr
->reppow
,12.0,10*GMX_DOUBLE_EPS
));
1283 fr
->bcoultab
= (!(fr
->eeltype
== eelCUT
|| EEL_RF(fr
->eeltype
)) ||
1284 fr
->eeltype
== eelRF_ZERO
);
1286 if (getenv("GMX_FORCE_TABLES"))
1289 fr
->bcoultab
= TRUE
;
1293 fprintf(fp
,"Table routines are used for coulomb: %s\n",bool_names
[fr
->bcoultab
]);
1294 fprintf(fp
,"Table routines are used for vdw: %s\n",bool_names
[fr
->bvdwtab
]);
1297 /* Tables are used for direct ewald sum */
1300 if (EEL_PME(ir
->coulombtype
))
1303 fprintf(fp
,"Will do PME sum in reciprocal space.\n");
1304 please_cite(fp
,"Essman95a");
1306 if (ir
->ewald_geometry
== eewg3DC
)
1310 fprintf(fp
,"Using the Ewald3DC correction for systems with a slab geometry.\n");
1312 please_cite(fp
,"In-Chul99a");
1315 fr
->ewaldcoeff
=calc_ewaldcoeff(ir
->rcoulomb
, ir
->ewald_rtol
);
1316 init_ewald_tab(&(fr
->ewald_table
), cr
, ir
, fp
);
1319 fprintf(fp
,"Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1324 /* Electrostatics */
1325 fr
->epsilon_r
= ir
->epsilon_r
;
1326 fr
->epsilon_rf
= ir
->epsilon_rf
;
1327 fr
->fudgeQQ
= mtop
->ffparams
.fudgeQQ
;
1328 fr
->rcoulomb_switch
= ir
->rcoulomb_switch
;
1329 fr
->rcoulomb
= cutoff_inf(ir
->rcoulomb
);
1331 /* Parameters for generalized RF */
1335 if (fr
->eeltype
== eelGRF
)
1337 init_generalized_rf(fp
,mtop
,ir
,fr
);
1339 else if (EEL_FULL(fr
->eeltype
) || (fr
->eeltype
== eelSHIFT
) ||
1340 (fr
->eeltype
== eelUSER
) || (fr
->eeltype
== eelSWITCH
))
1342 /* We must use the long range cut-off for neighboursearching...
1343 * An extra range of e.g. 0.1 nm (half the size of a charge group)
1344 * is necessary for neighboursearching. This allows diffusion
1345 * into the cut-off range (between neighborlist updates),
1346 * and gives more accurate forces because all atoms within the short-range
1347 * cut-off rc must be taken into account, while the ns criterium takes
1348 * only those with the center of geometry within the cut-off.
1349 * (therefore we have to add half the size of a charge group, plus
1350 * something to account for diffusion if we have nstlist > 1)
1352 for(m
=0; (m
<DIM
); m
++)
1353 box_size
[m
]=box
[m
][m
];
1355 if (fr
->eeltype
== eelPPPM
&& fr
->phi
== NULL
)
1356 snew(fr
->phi
,natoms
);
1358 if ((fr
->eeltype
==eelPPPM
) || (fr
->eeltype
==eelPOISSON
) ||
1359 (fr
->eeltype
== eelSHIFT
&& fr
->rcoulomb
> fr
->rcoulomb_switch
))
1360 set_shift_consts(fp
,fr
->rcoulomb_switch
,fr
->rcoulomb
,box_size
,fr
);
1363 fr
->bF_NoVirSum
= (EEL_FULL(fr
->eeltype
) ||
1364 gmx_mtop_ftype_count(mtop
,F_POSRES
) > 0 ||
1365 IR_ELEC_FIELD(*ir
));
1367 /* Mask that says whether or not this NBF list should be computed */
1368 /* if (fr->bMask == NULL) {
1369 ngrp = ir->opts.ngener*ir->opts.ngener;
1370 snew(fr->bMask,ngrp);*/
1371 /* Defaults to always */
1372 /* for(i=0; (i<ngrp); i++)
1373 fr->bMask[i] = TRUE;
1376 if (ncg_mtop(mtop
) > fr
->cg_nalloc
&& !DOMAINDECOMP(cr
)) {
1377 /* Count the total number of charge groups */
1378 fr
->cg_nalloc
= ncg_mtop(mtop
);
1379 srenew(fr
->cg_cm
,fr
->cg_nalloc
);
1381 if (fr
->shift_vec
== NULL
)
1382 snew(fr
->shift_vec
,SHIFTS
);
1384 if (fr
->fshift
== NULL
)
1385 snew(fr
->fshift
,SHIFTS
);
1387 if (fr
->nbfp
== NULL
) {
1388 fr
->ntype
= mtop
->ffparams
.atnr
;
1389 fr
->bBHAM
= (mtop
->ffparams
.functype
[0] == F_BHAM
);
1390 fr
->nbfp
= mk_nbfp(&mtop
->ffparams
,fr
->bBHAM
);
1393 /* Copy the energy group exclusions */
1394 fr
->egp_flags
= ir
->opts
.egp_flags
;
1396 /* Van der Waals stuff */
1397 fr
->rvdw
= cutoff_inf(ir
->rvdw
);
1398 fr
->rvdw_switch
= ir
->rvdw_switch
;
1399 if ((fr
->vdwtype
!= evdwCUT
) && (fr
->vdwtype
!= evdwUSER
) && !fr
->bBHAM
) {
1400 if (fr
->rvdw_switch
>= fr
->rvdw
)
1401 gmx_fatal(FARGS
,"rvdw_switch (%f) must be < rvdw (%f)",
1402 fr
->rvdw_switch
,fr
->rvdw
);
1404 fprintf(fp
,"Using %s Lennard-Jones, switch between %g and %g nm\n",
1405 (fr
->eeltype
==eelSWITCH
) ? "switched":"shifted",
1406 fr
->rvdw_switch
,fr
->rvdw
);
1409 if (fr
->bBHAM
&& (fr
->vdwtype
== evdwSHIFT
|| fr
->vdwtype
== evdwSWITCH
))
1410 gmx_fatal(FARGS
,"Switch/shift interaction not supported with Buckingham");
1413 fprintf(fp
,"Cut-off's: NS: %g Coulomb: %g %s: %g\n",
1414 fr
->rlist
,fr
->rcoulomb
,fr
->bBHAM
? "BHAM":"LJ",fr
->rvdw
);
1416 fr
->eDispCorr
= ir
->eDispCorr
;
1417 if (ir
->eDispCorr
!= edispcNO
)
1419 set_avcsixtwelve(fp
,fr
,mtop
);
1424 set_bham_b_max(fp
,fr
,mtop
);
1427 fr
->bGB
= (ir
->implicit_solvent
== eisGBSA
);
1428 fr
->gb_epsilon_solvent
= ir
->gb_epsilon_solvent
;
1430 /* Copy the GBSA data (radius, volume and surftens for each
1431 * atomtype) from the topology atomtype section to forcerec.
1433 snew(fr
->atype_radius
,fr
->ntype
);
1434 snew(fr
->atype_vol
,fr
->ntype
);
1435 snew(fr
->atype_surftens
,fr
->ntype
);
1436 snew(fr
->atype_gb_radius
,fr
->ntype
);
1437 snew(fr
->atype_S_hct
,fr
->ntype
);
1439 if (mtop
->atomtypes
.nr
> 0)
1441 for(i
=0;i
<fr
->ntype
;i
++)
1442 fr
->atype_radius
[i
] =mtop
->atomtypes
.radius
[i
];
1443 for(i
=0;i
<fr
->ntype
;i
++)
1444 fr
->atype_vol
[i
] = mtop
->atomtypes
.vol
[i
];
1445 for(i
=0;i
<fr
->ntype
;i
++)
1446 fr
->atype_surftens
[i
] = mtop
->atomtypes
.surftens
[i
];
1447 for(i
=0;i
<fr
->ntype
;i
++)
1448 fr
->atype_gb_radius
[i
] = mtop
->atomtypes
.gb_radius
[i
];
1449 for(i
=0;i
<fr
->ntype
;i
++)
1450 fr
->atype_S_hct
[i
] = mtop
->atomtypes
.S_hct
[i
];
1453 /* Generate the GB table if needed */
1457 fr
->gbtabscale
=2000;
1463 fr
->gbtab
=make_gb_table(fp
,oenv
,fr
,tabpfn
,fr
->gbtabscale
);
1465 init_gb(&fr
->born
,cr
,fr
,ir
,mtop
,ir
->rgbradii
,ir
->gb_algorithm
);
1467 /* Copy local gb data (for dd, this is done in dd_partition_system) */
1468 if (!DOMAINDECOMP(cr
))
1470 make_local_gb(cr
,fr
->born
,ir
->gb_algorithm
);
1474 /* Set the charge scaling */
1475 if (fr
->epsilon_r
!= 0)
1476 fr
->epsfac
= ONE_4PI_EPS0
/fr
->epsilon_r
;
1478 /* eps = 0 is infinite dieletric: no coulomb interactions */
1481 /* Reaction field constants */
1482 if (EEL_RF(fr
->eeltype
))
1483 calc_rffac(fp
,fr
->eeltype
,fr
->epsilon_r
,fr
->epsilon_rf
,
1484 fr
->rcoulomb
,fr
->temp
,fr
->zsquare
,box
,
1485 &fr
->kappa
,&fr
->k_rf
,&fr
->c_rf
);
1487 set_chargesum(fp
,fr
,mtop
);
1489 /* if we are using LR electrostatics, and they are tabulated,
1490 * the tables will contain modified coulomb interactions.
1491 * Since we want to use the non-shifted ones for 1-4
1492 * coulombic interactions, we must have an extra set of tables.
1495 /* Construct tables.
1496 * A little unnecessary to make both vdw and coul tables sometimes,
1497 * but what the heck... */
1499 bTab
= fr
->bcoultab
|| fr
->bvdwtab
;
1501 bSep14tab
= ((!bTab
|| fr
->eeltype
!=eelCUT
|| fr
->vdwtype
!=evdwCUT
||
1503 (gmx_mtop_ftype_count(mtop
,F_LJ14
) > 0 ||
1504 gmx_mtop_ftype_count(mtop
,F_LJC14_Q
) > 0 ||
1505 gmx_mtop_ftype_count(mtop
,F_LJC_PAIRS_NB
) > 0));
1507 negp_pp
= ir
->opts
.ngener
- ir
->nwall
;
1510 bNormalnblists
= TRUE
;
1513 bNormalnblists
= (ir
->eDispCorr
!= edispcNO
);
1514 for(egi
=0; egi
<negp_pp
; egi
++) {
1515 for(egj
=egi
; egj
<negp_pp
; egj
++) {
1516 egp_flags
= ir
->opts
.egp_flags
[GID(egi
,egj
,ir
->opts
.ngener
)];
1517 if (!(egp_flags
& EGP_EXCL
)) {
1518 if (egp_flags
& EGP_TABLE
) {
1521 bNormalnblists
= TRUE
;
1526 if (bNormalnblists
) {
1527 fr
->nnblists
= negptable
+ 1;
1529 fr
->nnblists
= negptable
;
1531 if (fr
->nnblists
> 1)
1532 snew(fr
->gid2nblists
,ir
->opts
.ngener
*ir
->opts
.ngener
);
1534 snew(fr
->nblists
,fr
->nnblists
);
1536 /* This code automatically gives table length tabext without cut-off's,
1537 * in that case grompp should already have checked that we do not need
1538 * normal tables and we only generate tables for 1-4 interactions.
1540 rtab
= ir
->rlistlong
+ ir
->tabext
;
1543 /* make tables for ordinary interactions */
1544 if (bNormalnblists
) {
1545 make_nbf_tables(fp
,oenv
,fr
,rtab
,cr
,tabfn
,NULL
,NULL
,&fr
->nblists
[0]);
1547 fr
->tab14
= fr
->nblists
[0].tab
;
1552 if (negptable
> 0) {
1553 /* Read the special tables for certain energy group pairs */
1554 nm_ind
= mtop
->groups
.grps
[egcENER
].nm_ind
;
1555 for(egi
=0; egi
<negp_pp
; egi
++) {
1556 for(egj
=egi
; egj
<negp_pp
; egj
++) {
1557 egp_flags
= ir
->opts
.egp_flags
[GID(egi
,egj
,ir
->opts
.ngener
)];
1558 if ((egp_flags
& EGP_TABLE
) && !(egp_flags
& EGP_EXCL
)) {
1559 nbl
= &(fr
->nblists
[m
]);
1560 if (fr
->nnblists
> 1) {
1561 fr
->gid2nblists
[GID(egi
,egj
,ir
->opts
.ngener
)] = m
;
1563 /* Read the table file with the two energy groups names appended */
1564 make_nbf_tables(fp
,oenv
,fr
,rtab
,cr
,tabfn
,
1565 *mtop
->groups
.grpname
[nm_ind
[egi
]],
1566 *mtop
->groups
.grpname
[nm_ind
[egj
]],
1569 } else if (fr
->nnblists
> 1) {
1570 fr
->gid2nblists
[GID(egi
,egj
,ir
->opts
.ngener
)] = 0;
1578 /* generate extra tables with plain Coulomb for 1-4 interactions only */
1579 fr
->tab14
= make_tables(fp
,oenv
,fr
,MASTER(cr
),tabpfn
,rtab
,
1580 GMX_MAKETABLES_14ONLY
);
1584 fr
->nwall
= ir
->nwall
;
1585 if (ir
->nwall
&& ir
->wall_type
==ewtTABLE
)
1587 make_wall_tables(fp
,oenv
,ir
,tabfn
,&mtop
->groups
,fr
);
1590 if (fcd
&& tabbfn
) {
1591 fcd
->bondtab
= make_bonded_tables(fp
,
1592 F_TABBONDS
,F_TABBONDSNC
,
1594 fcd
->angletab
= make_bonded_tables(fp
,
1597 fcd
->dihtab
= make_bonded_tables(fp
,
1602 fprintf(debug
,"No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
1605 if (ir
->eDispCorr
!= edispcNO
)
1607 calc_enervirdiff(fp
,ir
->eDispCorr
,fr
);
1610 /* QM/MM initialization if requested
1614 fprintf(stderr
,"QM/MM calculation requested.\n");
1617 fr
->bQMMM
= ir
->bQMMM
;
1618 fr
->qr
= mk_QMMMrec();
1620 /* Set all the static charge group info */
1621 fr
->cginfo_mb
= init_cginfo_mb(fp
,mtop
,fr
,bNoSolvOpt
);
1622 if (DOMAINDECOMP(cr
)) {
1625 fr
->cginfo
= cginfo_expand(mtop
->nmolblock
,fr
->cginfo_mb
);
1628 if (!DOMAINDECOMP(cr
))
1630 /* When using particle decomposition, the effect of the second argument,
1631 * which sets fr->hcg, is corrected later in do_md and init_em.
1633 forcerec_set_ranges(fr
,ncg_mtop(mtop
),ncg_mtop(mtop
),
1634 mtop
->natoms
,mtop
->natoms
);
1637 fr
->print_force
= print_force
;
1640 /* coarse load balancing vars */
1645 /* Initialize neighbor search */
1646 init_ns(fp
,cr
,&fr
->ns
,fr
,mtop
,box
);
1648 if (cr
->duty
& DUTY_PP
)
1649 gmx_setup_kernels(fp
);
1652 #define pr_real(fp,r) fprintf(fp,"%s: %e\n",#r,r)
1653 #define pr_int(fp,i) fprintf((fp),"%s: %d\n",#i,i)
1654 #define pr_bool(fp,b) fprintf((fp),"%s: %s\n",#b,bool_names[b])
1656 void pr_forcerec(FILE *fp
,t_forcerec
*fr
,t_commrec
*cr
)
1660 pr_real(fp
,fr
->rlist
);
1661 pr_real(fp
,fr
->rcoulomb
);
1662 pr_real(fp
,fr
->fudgeQQ
);
1663 pr_bool(fp
,fr
->bGrid
);
1664 pr_bool(fp
,fr
->bTwinRange
);
1665 /*pr_int(fp,fr->cg0);
1666 pr_int(fp,fr->hcg);*/
1667 for(i
=0; i
<fr
->nnblists
; i
++)
1668 pr_int(fp
,fr
->nblists
[i
].tab
.n
);
1669 pr_real(fp
,fr
->rcoulomb_switch
);
1670 pr_real(fp
,fr
->rcoulomb
);