2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/compat/make_unique.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald-utils.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/hw_info.h"
61 #include "gromacs/listed_forces/gpubonded.h"
62 #include "gromacs/listed_forces/manage-threading.h"
63 #include "gromacs/listed_forces/pairs.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/force.h"
69 #include "gromacs/mdlib/forcerec-threading.h"
70 #include "gromacs/mdlib/gmx_omp_nthreads.h"
71 #include "gromacs/mdlib/md_support.h"
72 #include "gromacs/mdlib/nb_verlet.h"
73 #include "gromacs/mdlib/nbnxn_atomdata.h"
74 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
75 #include "gromacs/mdlib/nbnxn_grid.h"
76 #include "gromacs/mdlib/nbnxn_internal.h"
77 #include "gromacs/mdlib/nbnxn_search.h"
78 #include "gromacs/mdlib/nbnxn_simd.h"
79 #include "gromacs/mdlib/nbnxn_tuning.h"
80 #include "gromacs/mdlib/nbnxn_util.h"
81 #include "gromacs/mdlib/ns.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/rf_util.h"
84 #include "gromacs/mdlib/sim_util.h"
85 #include "gromacs/mdlib/wall.h"
86 #include "gromacs/mdtypes/commrec.h"
87 #include "gromacs/mdtypes/fcdata.h"
88 #include "gromacs/mdtypes/group.h"
89 #include "gromacs/mdtypes/iforceprovider.h"
90 #include "gromacs/mdtypes/inputrec.h"
91 #include "gromacs/mdtypes/md_enums.h"
92 #include "gromacs/pbcutil/ishift.h"
93 #include "gromacs/pbcutil/pbc.h"
94 #include "gromacs/simd/simd.h"
95 #include "gromacs/tables/forcetable.h"
96 #include "gromacs/topology/mtop_util.h"
97 #include "gromacs/trajectory/trajectoryframe.h"
98 #include "gromacs/utility/cstringutil.h"
99 #include "gromacs/utility/exceptions.h"
100 #include "gromacs/utility/fatalerror.h"
101 #include "gromacs/utility/gmxassert.h"
102 #include "gromacs/utility/logger.h"
103 #include "gromacs/utility/physicalnodecommunicator.h"
104 #include "gromacs/utility/pleasecite.h"
105 #include "gromacs/utility/smalloc.h"
106 #include "gromacs/utility/strconvert.h"
108 #include "nbnxn_gpu_jit_support.h"
110 t_forcerec
*mk_forcerec()
119 static real
*mk_nbfp(const gmx_ffparams_t
*idef
, gmx_bool bBHAM
)
127 snew(nbfp
, 3*atnr
*atnr
);
128 for (i
= k
= 0; (i
< atnr
); i
++)
130 for (j
= 0; (j
< atnr
); j
++, k
++)
132 BHAMA(nbfp
, atnr
, i
, j
) = idef
->iparams
[k
].bham
.a
;
133 BHAMB(nbfp
, atnr
, i
, j
) = idef
->iparams
[k
].bham
.b
;
134 /* nbfp now includes the 6.0 derivative prefactor */
135 BHAMC(nbfp
, atnr
, i
, j
) = idef
->iparams
[k
].bham
.c
*6.0;
141 snew(nbfp
, 2*atnr
*atnr
);
142 for (i
= k
= 0; (i
< atnr
); i
++)
144 for (j
= 0; (j
< atnr
); j
++, k
++)
146 /* nbfp now includes the 6.0/12.0 derivative prefactors */
147 C6(nbfp
, atnr
, i
, j
) = idef
->iparams
[k
].lj
.c6
*6.0;
148 C12(nbfp
, atnr
, i
, j
) = idef
->iparams
[k
].lj
.c12
*12.0;
156 static real
*make_ljpme_c6grid(const gmx_ffparams_t
*idef
, t_forcerec
*fr
)
159 real c6
, c6i
, c6j
, c12i
, c12j
, epsi
, epsj
, sigmai
, sigmaj
;
162 /* For LJ-PME simulations, we correct the energies with the reciprocal space
163 * inside of the cut-off. To do this the non-bonded kernels needs to have
164 * access to the C6-values used on the reciprocal grid in pme.c
168 snew(grid
, 2*atnr
*atnr
);
169 for (i
= k
= 0; (i
< atnr
); i
++)
171 for (j
= 0; (j
< atnr
); j
++, k
++)
173 c6i
= idef
->iparams
[i
*(atnr
+1)].lj
.c6
;
174 c12i
= idef
->iparams
[i
*(atnr
+1)].lj
.c12
;
175 c6j
= idef
->iparams
[j
*(atnr
+1)].lj
.c6
;
176 c12j
= idef
->iparams
[j
*(atnr
+1)].lj
.c12
;
177 c6
= std::sqrt(c6i
* c6j
);
178 if (fr
->ljpme_combination_rule
== eljpmeLB
179 && !gmx_numzero(c6
) && !gmx_numzero(c12i
) && !gmx_numzero(c12j
))
181 sigmai
= gmx::sixthroot(c12i
/ c6i
);
182 sigmaj
= gmx::sixthroot(c12j
/ c6j
);
183 epsi
= c6i
* c6i
/ c12i
;
184 epsj
= c6j
* c6j
/ c12j
;
185 c6
= std::sqrt(epsi
* epsj
) * gmx::power6(0.5*(sigmai
+sigmaj
));
187 /* Store the elements at the same relative positions as C6 in nbfp in order
188 * to simplify access in the kernels
190 grid
[2*(atnr
*i
+j
)] = c6
*6.0;
196 static real
*mk_nbfp_combination_rule(const gmx_ffparams_t
*idef
, int comb_rule
)
200 real c6i
, c6j
, c12i
, c12j
, epsi
, epsj
, sigmai
, sigmaj
;
204 snew(nbfp
, 2*atnr
*atnr
);
205 for (i
= 0; i
< atnr
; ++i
)
207 for (j
= 0; j
< atnr
; ++j
)
209 c6i
= idef
->iparams
[i
*(atnr
+1)].lj
.c6
;
210 c12i
= idef
->iparams
[i
*(atnr
+1)].lj
.c12
;
211 c6j
= idef
->iparams
[j
*(atnr
+1)].lj
.c6
;
212 c12j
= idef
->iparams
[j
*(atnr
+1)].lj
.c12
;
213 c6
= std::sqrt(c6i
* c6j
);
214 c12
= std::sqrt(c12i
* c12j
);
215 if (comb_rule
== eCOMB_ARITHMETIC
216 && !gmx_numzero(c6
) && !gmx_numzero(c12
))
218 sigmai
= gmx::sixthroot(c12i
/ c6i
);
219 sigmaj
= gmx::sixthroot(c12j
/ c6j
);
220 epsi
= c6i
* c6i
/ c12i
;
221 epsj
= c6j
* c6j
/ c12j
;
222 c6
= std::sqrt(epsi
* epsj
) * gmx::power6(0.5*(sigmai
+sigmaj
));
223 c12
= std::sqrt(epsi
* epsj
) * gmx::power12(0.5*(sigmai
+sigmaj
));
225 C6(nbfp
, atnr
, i
, j
) = c6
*6.0;
226 C12(nbfp
, atnr
, i
, j
) = c12
*12.0;
232 /* This routine sets fr->solvent_opt to the most common solvent in the
233 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
234 * the fr->solvent_type array with the correct type (or esolNO).
236 * Charge groups that fulfill the conditions but are not identical to the
237 * most common one will be marked as esolNO in the solvent_type array.
239 * TIP3p is identical to SPC for these purposes, so we call it
240 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
242 * NOTE: QM particle should not
243 * become an optimized solvent. Not even if there is only one charge
253 } solvent_parameters_t
;
256 check_solvent_cg(const gmx_moltype_t
*molt
,
259 const unsigned char *qm_grpnr
,
260 const t_grps
*qm_grps
,
262 int *n_solvent_parameters
,
263 solvent_parameters_t
**solvent_parameters_p
,
273 real tmp_charge
[4] = { 0.0 }; /* init to zero to make gcc 7 happy */
274 int tmp_vdwtype
[4] = { 0 }; /* init to zero to make gcc 7 happy */
277 solvent_parameters_t
*solvent_parameters
;
279 /* We use a list with parameters for each solvent type.
280 * Every time we discover a new molecule that fulfills the basic
281 * conditions for a solvent we compare with the previous entries
282 * in these lists. If the parameters are the same we just increment
283 * the counter for that type, and otherwise we create a new type
284 * based on the current molecule.
286 * Once we've finished going through all molecules we check which
287 * solvent is most common, and mark all those molecules while we
288 * clear the flag on all others.
291 solvent_parameters
= *solvent_parameters_p
;
293 /* Mark the cg first as non optimized */
296 /* Check if this cg has no exclusions with atoms in other charge groups
297 * and all atoms inside the charge group excluded.
298 * We only have 3 or 4 atom solvent loops.
300 if (GET_CGINFO_EXCL_INTER(cginfo
) ||
301 !GET_CGINFO_EXCL_INTRA(cginfo
))
306 /* Get the indices of the first atom in this charge group */
307 j0
= molt
->cgs
.index
[cg0
];
308 j1
= molt
->cgs
.index
[cg0
+1];
310 /* Number of atoms in our molecule */
316 "Moltype '%s': there are %d atoms in this charge group\n",
320 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
323 if (nj
< 3 || nj
> 4)
328 /* Check if we are doing QM on this group */
330 if (qm_grpnr
!= nullptr)
332 for (j
= j0
; j
< j1
&& !qm
; j
++)
334 qm
= (qm_grpnr
[j
] < qm_grps
->nr
- 1);
337 /* Cannot use solvent optimization with QM */
343 atom
= molt
->atoms
.atom
;
345 /* Still looks like a solvent, time to check parameters */
347 /* If it is perturbed (free energy) we can't use the solvent loops,
348 * so then we just skip to the next molecule.
352 for (j
= j0
; j
< j1
&& !perturbed
; j
++)
354 perturbed
= PERTURBED(atom
[j
]);
362 /* Now it's only a question if the VdW and charge parameters
363 * are OK. Before doing the check we compare and see if they are
364 * identical to a possible previous solvent type.
365 * First we assign the current types and charges.
367 for (j
= 0; j
< nj
; j
++)
369 tmp_vdwtype
[j
] = atom
[j0
+j
].type
;
370 tmp_charge
[j
] = atom
[j0
+j
].q
;
373 /* Does it match any previous solvent type? */
374 for (k
= 0; k
< *n_solvent_parameters
; k
++)
379 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
380 if ( (solvent_parameters
[k
].model
== esolSPC
&& nj
!= 3) ||
381 (solvent_parameters
[k
].model
== esolTIP4P
&& nj
!= 4) )
386 /* Check that types & charges match for all atoms in molecule */
387 for (j
= 0; j
< nj
&& match
; j
++)
389 if (tmp_vdwtype
[j
] != solvent_parameters
[k
].vdwtype
[j
])
393 if (tmp_charge
[j
] != solvent_parameters
[k
].charge
[j
])
400 /* Congratulations! We have a matched solvent.
401 * Flag it with this type for later processing.
404 solvent_parameters
[k
].count
+= nmol
;
406 /* We are done with this charge group */
411 /* If we get here, we have a tentative new solvent type.
412 * Before we add it we must check that it fulfills the requirements
413 * of the solvent optimized loops. First determine which atoms have
416 for (j
= 0; j
< nj
; j
++)
419 tjA
= tmp_vdwtype
[j
];
421 /* Go through all other tpes and see if any have non-zero
422 * VdW parameters when combined with this one.
424 for (k
= 0; k
< fr
->ntype
&& (!has_vdw
[j
]); k
++)
426 /* We already checked that the atoms weren't perturbed,
427 * so we only need to check state A now.
431 has_vdw
[j
] = (has_vdw
[j
] ||
432 (BHAMA(fr
->nbfp
, fr
->ntype
, tjA
, k
) != 0.0) ||
433 (BHAMB(fr
->nbfp
, fr
->ntype
, tjA
, k
) != 0.0) ||
434 (BHAMC(fr
->nbfp
, fr
->ntype
, tjA
, k
) != 0.0));
439 has_vdw
[j
] = (has_vdw
[j
] ||
440 (C6(fr
->nbfp
, fr
->ntype
, tjA
, k
) != 0.0) ||
441 (C12(fr
->nbfp
, fr
->ntype
, tjA
, k
) != 0.0));
446 /* Now we know all we need to make the final check and assignment. */
450 * For this we require thatn all atoms have charge,
451 * the charges on atom 2 & 3 should be the same, and only
452 * atom 1 might have VdW.
456 tmp_charge
[0] != 0 &&
457 tmp_charge
[1] != 0 &&
458 tmp_charge
[2] == tmp_charge
[1])
460 srenew(solvent_parameters
, *n_solvent_parameters
+1);
461 solvent_parameters
[*n_solvent_parameters
].model
= esolSPC
;
462 solvent_parameters
[*n_solvent_parameters
].count
= nmol
;
463 for (k
= 0; k
< 3; k
++)
465 solvent_parameters
[*n_solvent_parameters
].vdwtype
[k
] = tmp_vdwtype
[k
];
466 solvent_parameters
[*n_solvent_parameters
].charge
[k
] = tmp_charge
[k
];
469 *cg_sp
= *n_solvent_parameters
;
470 (*n_solvent_parameters
)++;
475 /* Or could it be a TIP4P?
476 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
477 * Only atom 1 mght have VdW.
482 tmp_charge
[0] == 0 &&
483 tmp_charge
[1] != 0 &&
484 tmp_charge
[2] == tmp_charge
[1] &&
487 srenew(solvent_parameters
, *n_solvent_parameters
+1);
488 solvent_parameters
[*n_solvent_parameters
].model
= esolTIP4P
;
489 solvent_parameters
[*n_solvent_parameters
].count
= nmol
;
490 for (k
= 0; k
< 4; k
++)
492 solvent_parameters
[*n_solvent_parameters
].vdwtype
[k
] = tmp_vdwtype
[k
];
493 solvent_parameters
[*n_solvent_parameters
].charge
[k
] = tmp_charge
[k
];
496 *cg_sp
= *n_solvent_parameters
;
497 (*n_solvent_parameters
)++;
501 *solvent_parameters_p
= solvent_parameters
;
505 check_solvent(FILE * fp
,
506 const gmx_mtop_t
* mtop
,
508 cginfo_mb_t
*cginfo_mb
)
511 const gmx_moltype_t
*molt
;
512 int mol
, cg_mol
, at_offset
, am
, cgm
, i
, nmol_ch
, nmol
;
513 int n_solvent_parameters
;
514 solvent_parameters_t
*solvent_parameters
;
520 fprintf(debug
, "Going to determine what solvent types we have.\n");
523 n_solvent_parameters
= 0;
524 solvent_parameters
= nullptr;
525 /* Allocate temporary array for solvent type */
526 snew(cg_sp
, mtop
->molblock
.size());
529 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
531 molt
= &mtop
->moltype
[mtop
->molblock
[mb
].type
];
533 /* Here we have to loop over all individual molecules
534 * because we need to check for QMMM particles.
536 snew(cg_sp
[mb
], cginfo_mb
[mb
].cg_mod
);
537 nmol_ch
= cginfo_mb
[mb
].cg_mod
/cgs
->nr
;
538 nmol
= mtop
->molblock
[mb
].nmol
/nmol_ch
;
539 for (mol
= 0; mol
< nmol_ch
; mol
++)
542 am
= mol
*cgs
->index
[cgs
->nr
];
543 for (cg_mol
= 0; cg_mol
< cgs
->nr
; cg_mol
++)
545 check_solvent_cg(molt
, cg_mol
, nmol
,
546 mtop
->groups
.grpnr
[egcQMMM
] ?
547 mtop
->groups
.grpnr
[egcQMMM
]+at_offset
+am
: nullptr,
548 &mtop
->groups
.grps
[egcQMMM
],
550 &n_solvent_parameters
, &solvent_parameters
,
551 cginfo_mb
[mb
].cginfo
[cgm
+cg_mol
],
552 &cg_sp
[mb
][cgm
+cg_mol
]);
555 at_offset
+= cgs
->index
[cgs
->nr
];
558 /* Puh! We finished going through all charge groups.
559 * Now find the most common solvent model.
562 /* Most common solvent this far */
564 for (i
= 0; i
< n_solvent_parameters
; i
++)
567 solvent_parameters
[i
].count
> solvent_parameters
[bestsp
].count
)
575 bestsol
= solvent_parameters
[bestsp
].model
;
583 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
585 cgs
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].cgs
;
586 nmol
= (mtop
->molblock
[mb
].nmol
*cgs
->nr
)/cginfo_mb
[mb
].cg_mod
;
587 for (i
= 0; i
< cginfo_mb
[mb
].cg_mod
; i
++)
589 if (cg_sp
[mb
][i
] == bestsp
)
591 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[i
], bestsol
);
596 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[i
], esolNO
);
603 if (bestsol
!= esolNO
&& fp
!= nullptr)
605 fprintf(fp
, "\nEnabling %s-like water optimization for %d molecules.\n\n",
607 solvent_parameters
[bestsp
].count
);
610 sfree(solvent_parameters
);
611 fr
->solvent_opt
= bestsol
;
615 acNONE
= 0, acCONSTRAINT
, acSETTLE
618 static cginfo_mb_t
*init_cginfo_mb(FILE *fplog
, const gmx_mtop_t
*mtop
,
619 t_forcerec
*fr
, gmx_bool bNoSolvOpt
,
620 gmx_bool
*bFEP_NonBonded
,
621 gmx_bool
*bExcl_IntraCGAll_InterCGNone
)
624 const t_blocka
*excl
;
625 const gmx_moltype_t
*molt
;
626 const gmx_molblock_t
*molb
;
627 cginfo_mb_t
*cginfo_mb
;
630 int cg_offset
, a_offset
;
631 int m
, cg
, a0
, a1
, gid
, ai
, j
, aj
, excl_nalloc
;
635 gmx_bool bId
, *bExcl
, bExclIntraAll
, bExclInter
, bHaveVDW
, bHaveQ
, bHavePerturbedAtoms
;
637 snew(cginfo_mb
, mtop
->molblock
.size());
639 snew(type_VDW
, fr
->ntype
);
640 for (ai
= 0; ai
< fr
->ntype
; ai
++)
642 type_VDW
[ai
] = FALSE
;
643 for (j
= 0; j
< fr
->ntype
; j
++)
645 type_VDW
[ai
] = type_VDW
[ai
] ||
647 C6(fr
->nbfp
, fr
->ntype
, ai
, j
) != 0 ||
648 C12(fr
->nbfp
, fr
->ntype
, ai
, j
) != 0;
652 *bFEP_NonBonded
= FALSE
;
653 *bExcl_IntraCGAll_InterCGNone
= TRUE
;
656 snew(bExcl
, excl_nalloc
);
659 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
661 molb
= &mtop
->molblock
[mb
];
662 molt
= &mtop
->moltype
[molb
->type
];
666 /* Check if the cginfo is identical for all molecules in this block.
667 * If so, we only need an array of the size of one molecule.
668 * Otherwise we make an array of #mol times #cgs per molecule.
671 for (m
= 0; m
< molb
->nmol
; m
++)
673 int am
= m
*cgs
->index
[cgs
->nr
];
674 for (cg
= 0; cg
< cgs
->nr
; cg
++)
677 a1
= cgs
->index
[cg
+1];
678 if (getGroupType(mtop
->groups
, egcENER
, a_offset
+am
+a0
) !=
679 getGroupType(mtop
->groups
, egcENER
, a_offset
+a0
))
683 if (mtop
->groups
.grpnr
[egcQMMM
] != nullptr)
685 for (ai
= a0
; ai
< a1
; ai
++)
687 if (mtop
->groups
.grpnr
[egcQMMM
][a_offset
+am
+ai
] !=
688 mtop
->groups
.grpnr
[egcQMMM
][a_offset
+ai
])
697 cginfo_mb
[mb
].cg_start
= cg_offset
;
698 cginfo_mb
[mb
].cg_end
= cg_offset
+ molb
->nmol
*cgs
->nr
;
699 cginfo_mb
[mb
].cg_mod
= (bId
? 1 : molb
->nmol
)*cgs
->nr
;
700 snew(cginfo_mb
[mb
].cginfo
, cginfo_mb
[mb
].cg_mod
);
701 cginfo
= cginfo_mb
[mb
].cginfo
;
703 /* Set constraints flags for constrained atoms */
704 snew(a_con
, molt
->atoms
.nr
);
705 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
707 if (interaction_function
[ftype
].flags
& IF_CONSTRAINT
)
712 for (ia
= 0; ia
< molt
->ilist
[ftype
].size(); ia
+= 1+nral
)
716 for (a
= 0; a
< nral
; a
++)
718 a_con
[molt
->ilist
[ftype
].iatoms
[ia
+1+a
]] =
719 (ftype
== F_SETTLE
? acSETTLE
: acCONSTRAINT
);
725 for (m
= 0; m
< (bId
? 1 : molb
->nmol
); m
++)
728 int am
= m
*cgs
->index
[cgs
->nr
];
729 for (cg
= 0; cg
< cgs
->nr
; cg
++)
732 a1
= cgs
->index
[cg
+1];
734 /* Store the energy group in cginfo */
735 gid
= getGroupType(mtop
->groups
, egcENER
, a_offset
+am
+a0
);
736 SET_CGINFO_GID(cginfo
[cgm
+cg
], gid
);
738 /* Check the intra/inter charge group exclusions */
739 if (a1
-a0
> excl_nalloc
)
741 excl_nalloc
= a1
- a0
;
742 srenew(bExcl
, excl_nalloc
);
744 /* bExclIntraAll: all intra cg interactions excluded
745 * bExclInter: any inter cg interactions excluded
747 bExclIntraAll
= TRUE
;
751 bHavePerturbedAtoms
= FALSE
;
752 for (ai
= a0
; ai
< a1
; ai
++)
754 /* Check VDW and electrostatic interactions */
755 bHaveVDW
= bHaveVDW
|| (type_VDW
[molt
->atoms
.atom
[ai
].type
] ||
756 type_VDW
[molt
->atoms
.atom
[ai
].typeB
]);
757 bHaveQ
= bHaveQ
|| (molt
->atoms
.atom
[ai
].q
!= 0 ||
758 molt
->atoms
.atom
[ai
].qB
!= 0);
760 bHavePerturbedAtoms
= bHavePerturbedAtoms
|| (PERTURBED(molt
->atoms
.atom
[ai
]) != 0);
762 /* Clear the exclusion list for atom ai */
763 for (aj
= a0
; aj
< a1
; aj
++)
765 bExcl
[aj
-a0
] = FALSE
;
767 /* Loop over all the exclusions of atom ai */
768 for (j
= excl
->index
[ai
]; j
< excl
->index
[ai
+1]; j
++)
771 if (aj
< a0
|| aj
>= a1
)
780 /* Check if ai excludes a0 to a1 */
781 for (aj
= a0
; aj
< a1
; aj
++)
785 bExclIntraAll
= FALSE
;
792 SET_CGINFO_CONSTR(cginfo
[cgm
+cg
]);
795 SET_CGINFO_SETTLE(cginfo
[cgm
+cg
]);
803 SET_CGINFO_EXCL_INTRA(cginfo
[cgm
+cg
]);
807 SET_CGINFO_EXCL_INTER(cginfo
[cgm
+cg
]);
809 if (a1
- a0
> MAX_CHARGEGROUP_SIZE
)
811 /* The size in cginfo is currently only read with DD */
812 gmx_fatal(FARGS
, "A charge group has size %d which is larger than the limit of %d atoms", a1
-a0
, MAX_CHARGEGROUP_SIZE
);
816 SET_CGINFO_HAS_VDW(cginfo
[cgm
+cg
]);
820 SET_CGINFO_HAS_Q(cginfo
[cgm
+cg
]);
822 if (bHavePerturbedAtoms
&& fr
->efep
!= efepNO
)
824 SET_CGINFO_FEP(cginfo
[cgm
+cg
]);
825 *bFEP_NonBonded
= TRUE
;
827 /* Store the charge group size */
828 SET_CGINFO_NATOMS(cginfo
[cgm
+cg
], a1
-a0
);
830 if (!bExclIntraAll
|| bExclInter
)
832 *bExcl_IntraCGAll_InterCGNone
= FALSE
;
839 cg_offset
+= molb
->nmol
*cgs
->nr
;
840 a_offset
+= molb
->nmol
*cgs
->index
[cgs
->nr
];
845 /* the solvent optimizer is called after the QM is initialized,
846 * because we don't want to have the QM subsystemto become an
850 check_solvent(fplog
, mtop
, fr
, cginfo_mb
);
852 if (getenv("GMX_NO_SOLV_OPT"))
856 fprintf(fplog
, "Found environment variable GMX_NO_SOLV_OPT.\n"
857 "Disabling all solvent optimization\n");
859 fr
->solvent_opt
= esolNO
;
863 fr
->solvent_opt
= esolNO
;
865 if (!fr
->solvent_opt
)
867 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
869 for (cg
= 0; cg
< cginfo_mb
[mb
].cg_mod
; cg
++)
871 SET_CGINFO_SOLOPT(cginfo_mb
[mb
].cginfo
[cg
], esolNO
);
879 static int *cginfo_expand(int nmb
, cginfo_mb_t
*cgi_mb
)
884 ncg
= cgi_mb
[nmb
-1].cg_end
;
887 for (cg
= 0; cg
< ncg
; cg
++)
889 while (cg
>= cgi_mb
[mb
].cg_end
)
894 cgi_mb
[mb
].cginfo
[(cg
- cgi_mb
[mb
].cg_start
) % cgi_mb
[mb
].cg_mod
];
900 static void done_cginfo_mb(cginfo_mb_t
*cginfo_mb
, int numMolBlocks
)
902 if (cginfo_mb
== nullptr)
906 for (int mb
= 0; mb
< numMolBlocks
; ++mb
)
908 sfree(cginfo_mb
[mb
].cginfo
);
913 /* Sets the sum of charges (squared) and C6 in the system in fr.
914 * Returns whether the system has a net charge.
916 static bool set_chargesum(FILE *log
, t_forcerec
*fr
, const gmx_mtop_t
*mtop
)
918 /*This now calculates sum for q and c6*/
919 double qsum
, q2sum
, q
, c6sum
, c6
;
924 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
926 int nmol
= molb
.nmol
;
927 const t_atoms
*atoms
= &mtop
->moltype
[molb
.type
].atoms
;
928 for (int i
= 0; i
< atoms
->nr
; i
++)
930 q
= atoms
->atom
[i
].q
;
933 c6
= mtop
->ffparams
.iparams
[atoms
->atom
[i
].type
*(mtop
->ffparams
.atnr
+1)].lj
.c6
;
938 fr
->q2sum
[0] = q2sum
;
939 fr
->c6sum
[0] = c6sum
;
941 if (fr
->efep
!= efepNO
)
946 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
948 int nmol
= molb
.nmol
;
949 const t_atoms
*atoms
= &mtop
->moltype
[molb
.type
].atoms
;
950 for (int i
= 0; i
< atoms
->nr
; i
++)
952 q
= atoms
->atom
[i
].qB
;
955 c6
= mtop
->ffparams
.iparams
[atoms
->atom
[i
].typeB
*(mtop
->ffparams
.atnr
+1)].lj
.c6
;
959 fr
->q2sum
[1] = q2sum
;
960 fr
->c6sum
[1] = c6sum
;
965 fr
->qsum
[1] = fr
->qsum
[0];
966 fr
->q2sum
[1] = fr
->q2sum
[0];
967 fr
->c6sum
[1] = fr
->c6sum
[0];
971 if (fr
->efep
== efepNO
)
973 fprintf(log
, "System total charge: %.3f\n", fr
->qsum
[0]);
977 fprintf(log
, "System total charge, top. A: %.3f top. B: %.3f\n",
978 fr
->qsum
[0], fr
->qsum
[1]);
982 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
983 return (std::abs(fr
->qsum
[0]) > 1e-4 ||
984 std::abs(fr
->qsum
[1]) > 1e-4);
987 void update_forcerec(t_forcerec
*fr
, matrix box
)
989 if (fr
->ic
->eeltype
== eelGRF
)
991 calc_rffac(nullptr, fr
->ic
->eeltype
, fr
->ic
->epsilon_r
, fr
->ic
->epsilon_rf
,
992 fr
->ic
->rcoulomb
, fr
->temp
, fr
->zsquare
, box
,
993 &fr
->ic
->k_rf
, &fr
->ic
->c_rf
);
997 void set_avcsixtwelve(FILE *fplog
, t_forcerec
*fr
, const gmx_mtop_t
*mtop
)
999 const t_atoms
*atoms
, *atoms_tpi
;
1000 const t_blocka
*excl
;
1001 int nmolc
, i
, j
, tpi
, tpj
, j1
, j2
, k
, nexcl
, q
;
1002 int64_t npair
, npair_ij
, tmpi
, tmpj
;
1003 double csix
, ctwelve
;
1004 int ntp
, *typecount
;
1007 real
*nbfp_comb
= nullptr;
1013 /* For LJ-PME, we want to correct for the difference between the
1014 * actual C6 values and the C6 values used by the LJ-PME based on
1015 * combination rules. */
1017 if (EVDW_PME(fr
->ic
->vdwtype
))
1019 nbfp_comb
= mk_nbfp_combination_rule(&mtop
->ffparams
,
1020 (fr
->ljpme_combination_rule
== eljpmeLB
) ? eCOMB_ARITHMETIC
: eCOMB_GEOMETRIC
);
1021 for (tpi
= 0; tpi
< ntp
; ++tpi
)
1023 for (tpj
= 0; tpj
< ntp
; ++tpj
)
1025 C6(nbfp_comb
, ntp
, tpi
, tpj
) =
1026 C6(nbfp
, ntp
, tpi
, tpj
) - C6(nbfp_comb
, ntp
, tpi
, tpj
);
1027 C12(nbfp_comb
, ntp
, tpi
, tpj
) = C12(nbfp
, ntp
, tpi
, tpj
);
1032 for (q
= 0; q
< (fr
->efep
== efepNO
? 1 : 2); q
++)
1040 /* Count the types so we avoid natoms^2 operations */
1041 snew(typecount
, ntp
);
1042 gmx_mtop_count_atomtypes(mtop
, q
, typecount
);
1044 for (tpi
= 0; tpi
< ntp
; tpi
++)
1046 for (tpj
= tpi
; tpj
< ntp
; tpj
++)
1048 tmpi
= typecount
[tpi
];
1049 tmpj
= typecount
[tpj
];
1052 npair_ij
= tmpi
*tmpj
;
1056 npair_ij
= tmpi
*(tmpi
- 1)/2;
1060 /* nbfp now includes the 6.0 derivative prefactor */
1061 csix
+= npair_ij
*BHAMC(nbfp
, ntp
, tpi
, tpj
)/6.0;
1065 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1066 csix
+= npair_ij
* C6(nbfp
, ntp
, tpi
, tpj
)/6.0;
1067 ctwelve
+= npair_ij
* C12(nbfp
, ntp
, tpi
, tpj
)/12.0;
1073 /* Subtract the excluded pairs.
1074 * The main reason for substracting exclusions is that in some cases
1075 * some combinations might never occur and the parameters could have
1076 * any value. These unused values should not influence the dispersion
1079 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
1081 int nmol
= molb
.nmol
;
1082 atoms
= &mtop
->moltype
[molb
.type
].atoms
;
1083 excl
= &mtop
->moltype
[molb
.type
].excls
;
1084 for (int i
= 0; (i
< atoms
->nr
); i
++)
1088 tpi
= atoms
->atom
[i
].type
;
1092 tpi
= atoms
->atom
[i
].typeB
;
1094 j1
= excl
->index
[i
];
1095 j2
= excl
->index
[i
+1];
1096 for (j
= j1
; j
< j2
; j
++)
1103 tpj
= atoms
->atom
[k
].type
;
1107 tpj
= atoms
->atom
[k
].typeB
;
1111 /* nbfp now includes the 6.0 derivative prefactor */
1112 csix
-= nmol
*BHAMC(nbfp
, ntp
, tpi
, tpj
)/6.0;
1116 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1117 csix
-= nmol
*C6 (nbfp
, ntp
, tpi
, tpj
)/6.0;
1118 ctwelve
-= nmol
*C12(nbfp
, ntp
, tpi
, tpj
)/12.0;
1128 /* Only correct for the interaction of the test particle
1129 * with the rest of the system.
1132 &mtop
->moltype
[mtop
->molblock
.back().type
].atoms
;
1135 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
1137 const gmx_molblock_t
&molb
= mtop
->molblock
[mb
];
1138 atoms
= &mtop
->moltype
[molb
.type
].atoms
;
1139 for (j
= 0; j
< atoms
->nr
; j
++)
1142 /* Remove the interaction of the test charge group
1145 if (mb
== mtop
->molblock
.size() - 1)
1149 if (mb
== 0 && molb
.nmol
== 1)
1151 gmx_fatal(FARGS
, "Old format tpr with TPI, please generate a new tpr file");
1156 tpj
= atoms
->atom
[j
].type
;
1160 tpj
= atoms
->atom
[j
].typeB
;
1162 for (i
= 0; i
< fr
->n_tpi
; i
++)
1166 tpi
= atoms_tpi
->atom
[i
].type
;
1170 tpi
= atoms_tpi
->atom
[i
].typeB
;
1174 /* nbfp now includes the 6.0 derivative prefactor */
1175 csix
+= nmolc
*BHAMC(nbfp
, ntp
, tpi
, tpj
)/6.0;
1179 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1180 csix
+= nmolc
*C6 (nbfp
, ntp
, tpi
, tpj
)/6.0;
1181 ctwelve
+= nmolc
*C12(nbfp
, ntp
, tpi
, tpj
)/12.0;
1188 if (npair
- nexcl
<= 0 && fplog
)
1190 fprintf(fplog
, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1196 csix
/= npair
- nexcl
;
1197 ctwelve
/= npair
- nexcl
;
1201 fprintf(debug
, "Counted %d exclusions\n", nexcl
);
1202 fprintf(debug
, "Average C6 parameter is: %10g\n", csix
);
1203 fprintf(debug
, "Average C12 parameter is: %10g\n", ctwelve
);
1205 fr
->avcsix
[q
] = csix
;
1206 fr
->avctwelve
[q
] = ctwelve
;
1209 if (EVDW_PME(fr
->ic
->vdwtype
))
1214 if (fplog
!= nullptr)
1216 if (fr
->eDispCorr
== edispcAllEner
||
1217 fr
->eDispCorr
== edispcAllEnerPres
)
1219 fprintf(fplog
, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1220 fr
->avcsix
[0], fr
->avctwelve
[0]);
1224 fprintf(fplog
, "Long Range LJ corr.: <C6> %10.4e\n", fr
->avcsix
[0]);
1230 static real
calcBuckinghamBMax(FILE *fplog
, const gmx_mtop_t
*mtop
)
1232 const t_atoms
*at1
, *at2
;
1233 int i
, j
, tpi
, tpj
, ntypes
;
1238 fprintf(fplog
, "Determining largest Buckingham b parameter for table\n");
1240 ntypes
= mtop
->ffparams
.atnr
;
1243 real bham_b_max
= 0;
1244 for (size_t mt1
= 0; mt1
< mtop
->moltype
.size(); mt1
++)
1246 at1
= &mtop
->moltype
[mt1
].atoms
;
1247 for (i
= 0; (i
< at1
->nr
); i
++)
1249 tpi
= at1
->atom
[i
].type
;
1252 gmx_fatal(FARGS
, "Atomtype[%d] = %d, maximum = %d", i
, tpi
, ntypes
);
1255 for (size_t mt2
= mt1
; mt2
< mtop
->moltype
.size(); mt2
++)
1257 at2
= &mtop
->moltype
[mt2
].atoms
;
1258 for (j
= 0; (j
< at2
->nr
); j
++)
1260 tpj
= at2
->atom
[j
].type
;
1263 gmx_fatal(FARGS
, "Atomtype[%d] = %d, maximum = %d", j
, tpj
, ntypes
);
1265 b
= mtop
->ffparams
.iparams
[tpi
*ntypes
+ tpj
].bham
.b
;
1270 if ((b
< bmin
) || (bmin
== -1))
1280 fprintf(fplog
, "Buckingham b parameters, min: %g, max: %g\n",
1287 static void make_nbf_tables(FILE *fp
,
1288 const interaction_const_t
*ic
, real rtab
,
1289 const char *tabfn
, char *eg1
, char *eg2
,
1295 if (tabfn
== nullptr)
1299 fprintf(debug
, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1304 sprintf(buf
, "%s", tabfn
);
1307 /* Append the two energy group names */
1308 sprintf(buf
+ strlen(tabfn
) - strlen(ftp2ext(efXVG
)) - 1, "_%s_%s.%s",
1309 eg1
, eg2
, ftp2ext(efXVG
));
1311 nbl
->table_elec_vdw
= make_tables(fp
, ic
, buf
, rtab
, 0);
1312 /* Copy the contents of the table to separate coulomb and LJ tables too,
1313 * to improve cache performance.
1315 /* For performance reasons we want
1316 * the table data to be aligned to 16-byte. The pointers could be freed
1317 * but currently aren't.
1319 snew(nbl
->table_elec
, 1);
1320 nbl
->table_elec
->interaction
= GMX_TABLE_INTERACTION_ELEC
;
1321 nbl
->table_elec
->format
= nbl
->table_elec_vdw
->format
;
1322 nbl
->table_elec
->r
= nbl
->table_elec_vdw
->r
;
1323 nbl
->table_elec
->n
= nbl
->table_elec_vdw
->n
;
1324 nbl
->table_elec
->scale
= nbl
->table_elec_vdw
->scale
;
1325 nbl
->table_elec
->formatsize
= nbl
->table_elec_vdw
->formatsize
;
1326 nbl
->table_elec
->ninteractions
= 1;
1327 nbl
->table_elec
->stride
= nbl
->table_elec
->formatsize
* nbl
->table_elec
->ninteractions
;
1328 snew_aligned(nbl
->table_elec
->data
, nbl
->table_elec
->stride
*(nbl
->table_elec
->n
+1), 32);
1330 snew(nbl
->table_vdw
, 1);
1331 nbl
->table_vdw
->interaction
= GMX_TABLE_INTERACTION_VDWREP_VDWDISP
;
1332 nbl
->table_vdw
->format
= nbl
->table_elec_vdw
->format
;
1333 nbl
->table_vdw
->r
= nbl
->table_elec_vdw
->r
;
1334 nbl
->table_vdw
->n
= nbl
->table_elec_vdw
->n
;
1335 nbl
->table_vdw
->scale
= nbl
->table_elec_vdw
->scale
;
1336 nbl
->table_vdw
->formatsize
= nbl
->table_elec_vdw
->formatsize
;
1337 nbl
->table_vdw
->ninteractions
= 2;
1338 nbl
->table_vdw
->stride
= nbl
->table_vdw
->formatsize
* nbl
->table_vdw
->ninteractions
;
1339 snew_aligned(nbl
->table_vdw
->data
, nbl
->table_vdw
->stride
*(nbl
->table_vdw
->n
+1), 32);
1341 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1342 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1344 for (i
= 0; i
<= nbl
->table_elec_vdw
->n
; i
++)
1346 for (j
= 0; j
< 4; j
++)
1348 nbl
->table_elec
->data
[4*i
+j
] = nbl
->table_elec_vdw
->data
[12*i
+j
];
1351 for (i
= 0; i
<= nbl
->table_elec_vdw
->n
; i
++)
1353 for (j
= 0; j
< 8; j
++)
1355 nbl
->table_vdw
->data
[8*i
+j
] = nbl
->table_elec_vdw
->data
[12*i
+4+j
];
1360 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1361 * ftype2 present in the topology, build an array of the number of
1362 * interactions present for each bonded interaction index found in the
1365 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1366 * valid type with that parameter.
1368 * \c count will be reallocated as necessary to fit the largest bonded
1369 * interaction index found, and its current size will be returned in
1370 * \c ncount. It will contain zero for every bonded interaction index
1371 * for which no interactions are present in the topology.
1373 static void count_tables(int ftype1
, int ftype2
, const gmx_mtop_t
*mtop
,
1374 int *ncount
, int **count
)
1376 int ftype
, i
, j
, tabnr
;
1378 // Loop over all moleculetypes
1379 for (const gmx_moltype_t
&molt
: mtop
->moltype
)
1381 // Loop over all interaction types
1382 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1384 // If the current interaction type is one of the types whose tables we're trying to count...
1385 if (ftype
== ftype1
|| ftype
== ftype2
)
1387 const InteractionList
&il
= molt
.ilist
[ftype
];
1388 const int stride
= 1 + NRAL(ftype
);
1389 // ... and there are actually some interactions for this type
1390 for (i
= 0; i
< il
.size(); i
+= stride
)
1392 // Find out which table index the user wanted
1393 tabnr
= mtop
->ffparams
.iparams
[il
.iatoms
[i
]].tab
.table
;
1396 gmx_fatal(FARGS
, "A bonded table number is smaller than 0: %d\n", tabnr
);
1398 // Make room for this index in the data structure
1399 if (tabnr
>= *ncount
)
1401 srenew(*count
, tabnr
+1);
1402 for (j
= *ncount
; j
< tabnr
+1; j
++)
1408 // Record that this table index is used and must have a valid file
1416 /*!\brief If there's bonded interactions of flavour \c tabext and type
1417 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1418 * list of filenames passed to mdrun, and make bonded tables from
1421 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1422 * valid type with that parameter.
1424 * A fatal error occurs if no matching filename is found.
1426 static bondedtable_t
*make_bonded_tables(FILE *fplog
,
1427 int ftype1
, int ftype2
,
1428 const gmx_mtop_t
*mtop
,
1429 gmx::ArrayRef
<const std::string
> tabbfnm
,
1439 count_tables(ftype1
, ftype2
, mtop
, &ncount
, &count
);
1441 // Are there any relevant tabulated bond interactions?
1445 for (int i
= 0; i
< ncount
; i
++)
1447 // Do any interactions exist that requires this table?
1450 // This pattern enforces the current requirement that
1451 // table filenames end in a characteristic sequence
1452 // before the file type extension, and avoids table 13
1453 // being recognized and used for table 1.
1454 std::string patternToFind
= gmx::formatString("_%s%d.%s", tabext
, i
, ftp2ext(efXVG
));
1455 bool madeTable
= false;
1456 for (gmx::index j
= 0; j
< tabbfnm
.size() && !madeTable
; ++j
)
1458 if (gmx::endsWith(tabbfnm
[j
], patternToFind
))
1460 // Finally read the table from the file found
1461 tab
[i
] = make_bonded_table(fplog
, tabbfnm
[j
].c_str(), NRAL(ftype1
)-2);
1467 bool isPlural
= (ftype2
!= -1);
1468 gmx_fatal(FARGS
, "Tabulated interaction of type '%s%s%s' with index %d cannot be used because no table file whose name matched '%s' was passed via the gmx mdrun -tableb command-line option.",
1469 interaction_function
[ftype1
].longname
,
1470 isPlural
? "' or '" : "",
1471 isPlural
? interaction_function
[ftype2
].longname
: "",
1473 patternToFind
.c_str());
1483 void forcerec_set_ranges(t_forcerec
*fr
,
1484 int ncg_home
, int ncg_force
,
1486 int natoms_force_constr
, int natoms_f_novirsum
)
1491 /* fr->ncg_force is unused in the standard code,
1492 * but it can be useful for modified code dealing with charge groups.
1494 fr
->ncg_force
= ncg_force
;
1495 fr
->natoms_force
= natoms_force
;
1496 fr
->natoms_force_constr
= natoms_force_constr
;
1498 if (fr
->natoms_force_constr
> fr
->nalloc_force
)
1500 fr
->nalloc_force
= over_alloc_dd(fr
->natoms_force_constr
);
1503 if (fr
->haveDirectVirialContributions
)
1505 fr
->forceBufferForDirectVirialContributions
->resize(natoms_f_novirsum
);
1509 static real
cutoff_inf(real cutoff
)
1513 cutoff
= GMX_CUTOFF_INF
;
1519 gmx_bool
nbnxn_simd_supported(const gmx::MDLogger
&mdlog
,
1520 const t_inputrec
*ir
)
1522 if (ir
->vdwtype
== evdwPME
&& ir
->ljpme_combination_rule
== eljpmeLB
)
1524 /* LJ PME with LB combination rule does 7 mesh operations.
1525 * This so slow that we don't compile SIMD non-bonded kernels
1527 GMX_LOG(mdlog
.warning
).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1535 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused
*ir
,
1538 const gmx_hw_info_t gmx_unused
&hardwareInfo
)
1540 *kernel_type
= nbnxnk4x4_PlainC
;
1541 *ewald_excl
= ewaldexclTable
;
1545 #ifdef GMX_NBNXN_SIMD_4XN
1546 *kernel_type
= nbnxnk4xN_SIMD_4xN
;
1548 #ifdef GMX_NBNXN_SIMD_2XNN
1549 *kernel_type
= nbnxnk4xN_SIMD_2xNN
;
1552 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1553 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1554 * This is based on the SIMD acceleration choice and CPU information
1555 * detected at runtime.
1557 * 4xN calculates more (zero) interactions, but has less pair-search
1558 * work and much better kernel instruction scheduling.
1560 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1561 * which doesn't have FMA, both the analytical and tabulated Ewald
1562 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1563 * 2x(4+4) because it results in significantly fewer pairs.
1564 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1565 * 10% with HT, 50% without HT. As we currently don't detect the actual
1566 * use of HT, use 4x8 to avoid a potential performance hit.
1567 * On Intel Haswell 4x8 is always faster.
1569 *kernel_type
= nbnxnk4xN_SIMD_4xN
;
1571 #if !GMX_SIMD_HAVE_FMA
1572 if (EEL_PME_EWALD(ir
->coulombtype
) ||
1573 EVDW_PME(ir
->vdwtype
))
1575 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1576 * There are enough instructions to make 2x(4+4) efficient.
1578 *kernel_type
= nbnxnk4xN_SIMD_2xNN
;
1581 if (hardwareInfo
.haveAmdZenCpu
)
1583 /* One 256-bit FMA per cycle makes 2xNN faster */
1584 *kernel_type
= nbnxnk4xN_SIMD_2xNN
;
1586 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1589 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1591 #ifdef GMX_NBNXN_SIMD_4XN
1592 *kernel_type
= nbnxnk4xN_SIMD_4xN
;
1594 gmx_fatal(FARGS
, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1597 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1599 #ifdef GMX_NBNXN_SIMD_2XNN
1600 *kernel_type
= nbnxnk4xN_SIMD_2xNN
;
1602 gmx_fatal(FARGS
, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1606 /* Analytical Ewald exclusion correction is only an option in
1608 * Since table lookup's don't parallelize with SIMD, analytical
1609 * will probably always be faster for a SIMD width of 8 or more.
1610 * With FMA analytical is sometimes faster for a width if 4 as well.
1611 * In single precision, this is faster on Bulldozer.
1613 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1614 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1615 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1616 * of single or double precision and 128 or 256-bit AVX2.
1618 if (!hardwareInfo
.haveAmdZenCpu
)
1620 *ewald_excl
= ewaldexclAnalytical
;
1623 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1625 *ewald_excl
= ewaldexclTable
;
1627 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1629 *ewald_excl
= ewaldexclAnalytical
;
1637 const char *lookup_nbnxn_kernel_name(int kernel_type
)
1639 const char *returnvalue
= nullptr;
1640 switch (kernel_type
)
1643 returnvalue
= "not set";
1645 case nbnxnk4x4_PlainC
:
1646 returnvalue
= "plain C";
1648 case nbnxnk4xN_SIMD_4xN
:
1649 case nbnxnk4xN_SIMD_2xNN
:
1651 returnvalue
= "SIMD";
1653 returnvalue
= "not available";
1656 case nbnxnk8x8x8_GPU
: returnvalue
= "GPU"; break;
1657 case nbnxnk8x8x8_PlainC
: returnvalue
= "plain C"; break;
1661 gmx_fatal(FARGS
, "Illegal kernel type selected");
1666 static void pick_nbnxn_kernel(const gmx::MDLogger
&mdlog
,
1667 gmx_bool use_simd_kernels
,
1668 const gmx_hw_info_t
&hardwareInfo
,
1670 EmulateGpuNonbonded emulateGpu
,
1671 const t_inputrec
*ir
,
1674 gmx_bool bDoNonbonded
)
1676 assert(kernel_type
);
1678 *kernel_type
= nbnxnkNotSet
;
1679 *ewald_excl
= ewaldexclTable
;
1681 if (emulateGpu
== EmulateGpuNonbonded::Yes
)
1683 *kernel_type
= nbnxnk8x8x8_PlainC
;
1687 GMX_LOG(mdlog
.warning
).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1692 *kernel_type
= nbnxnk8x8x8_GPU
;
1695 if (*kernel_type
== nbnxnkNotSet
)
1697 if (use_simd_kernels
&&
1698 nbnxn_simd_supported(mdlog
, ir
))
1700 pick_nbnxn_kernel_cpu(ir
, kernel_type
, ewald_excl
, hardwareInfo
);
1704 *kernel_type
= nbnxnk4x4_PlainC
;
1710 GMX_LOG(mdlog
.info
).asParagraph().appendTextFormatted(
1711 "Using %s %dx%d nonbonded short-range kernels",
1712 lookup_nbnxn_kernel_name(*kernel_type
),
1713 nbnxn_kernel_to_cluster_i_size(*kernel_type
),
1714 nbnxn_kernel_to_cluster_j_size(*kernel_type
));
1716 if (nbnxnk4x4_PlainC
== *kernel_type
||
1717 nbnxnk8x8x8_PlainC
== *kernel_type
)
1719 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1720 "WARNING: Using the slow %s kernels. This should\n"
1721 "not happen during routine usage on supported platforms.",
1722 lookup_nbnxn_kernel_name(*kernel_type
));
1727 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1728 static void initCoulombEwaldParameters(FILE *fp
, const t_inputrec
*ir
,
1729 bool systemHasNetCharge
,
1730 interaction_const_t
*ic
)
1732 if (!EEL_PME_EWALD(ir
->coulombtype
))
1739 fprintf(fp
, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1741 if (ir
->coulombtype
== eelP3M_AD
)
1743 please_cite(fp
, "Hockney1988");
1744 please_cite(fp
, "Ballenegger2012");
1748 please_cite(fp
, "Essmann95a");
1751 if (ir
->ewald_geometry
== eewg3DC
)
1755 fprintf(fp
, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1756 systemHasNetCharge
? " and net charge" : "");
1758 please_cite(fp
, "In-Chul99a");
1759 if (systemHasNetCharge
)
1761 please_cite(fp
, "Ballenegger2009");
1766 ic
->ewaldcoeff_q
= calc_ewaldcoeff_q(ir
->rcoulomb
, ir
->ewald_rtol
);
1769 fprintf(fp
, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1770 1/ic
->ewaldcoeff_q
);
1773 if (ic
->coulomb_modifier
== eintmodPOTSHIFT
)
1775 GMX_RELEASE_ASSERT(ic
->rcoulomb
!= 0, "Cutoff radius cannot be zero");
1776 ic
->sh_ewald
= std::erfc(ic
->ewaldcoeff_q
*ic
->rcoulomb
) / ic
->rcoulomb
;
1784 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1785 static void initVdwEwaldParameters(FILE *fp
, const t_inputrec
*ir
,
1786 interaction_const_t
*ic
)
1788 if (!EVDW_PME(ir
->vdwtype
))
1795 fprintf(fp
, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1796 please_cite(fp
, "Essmann95a");
1798 ic
->ewaldcoeff_lj
= calc_ewaldcoeff_lj(ir
->rvdw
, ir
->ewald_rtol_lj
);
1801 fprintf(fp
, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1802 1/ic
->ewaldcoeff_lj
);
1805 if (ic
->vdw_modifier
== eintmodPOTSHIFT
)
1807 real crc2
= gmx::square(ic
->ewaldcoeff_lj
*ic
->rvdw
);
1808 ic
->sh_lj_ewald
= (std::exp(-crc2
)*(1 + crc2
+ 0.5*crc2
*crc2
) - 1)/gmx::power6(ic
->rvdw
);
1812 ic
->sh_lj_ewald
= 0;
1816 gmx_bool
uses_simple_tables(int cutoff_scheme
,
1817 nonbonded_verlet_t
*nbv
,
1820 gmx_bool bUsesSimpleTables
= TRUE
;
1823 switch (cutoff_scheme
)
1826 bUsesSimpleTables
= TRUE
;
1829 assert(nullptr != nbv
);
1830 grp_index
= (group
< 0) ? 0 : (nbv
->ngrp
- 1);
1831 bUsesSimpleTables
= nbnxn_kernel_pairlist_simple(nbv
->grp
[grp_index
].kernel_type
);
1834 gmx_incons("unimplemented");
1836 return bUsesSimpleTables
;
1839 static void init_ewald_f_table(interaction_const_t
*ic
,
1844 /* Get the Ewald table spacing based on Coulomb and/or LJ
1845 * Ewald coefficients and rtol.
1847 ic
->tabq_scale
= ewald_spline3_table_scale(ic
);
1849 if (ic
->cutoff_scheme
== ecutsVERLET
)
1851 maxr
= ic
->rcoulomb
;
1855 maxr
= std::max(ic
->rcoulomb
, rtab
);
1857 ic
->tabq_size
= static_cast<int>(maxr
*ic
->tabq_scale
) + 2;
1859 sfree_aligned(ic
->tabq_coul_FDV0
);
1860 sfree_aligned(ic
->tabq_coul_F
);
1861 sfree_aligned(ic
->tabq_coul_V
);
1863 sfree_aligned(ic
->tabq_vdw_FDV0
);
1864 sfree_aligned(ic
->tabq_vdw_F
);
1865 sfree_aligned(ic
->tabq_vdw_V
);
1867 if (EEL_PME_EWALD(ic
->eeltype
))
1869 /* Create the original table data in FDV0 */
1870 snew_aligned(ic
->tabq_coul_FDV0
, ic
->tabq_size
*4, 32);
1871 snew_aligned(ic
->tabq_coul_F
, ic
->tabq_size
, 32);
1872 snew_aligned(ic
->tabq_coul_V
, ic
->tabq_size
, 32);
1873 table_spline3_fill_ewald_lr(ic
->tabq_coul_F
, ic
->tabq_coul_V
, ic
->tabq_coul_FDV0
,
1874 ic
->tabq_size
, 1/ic
->tabq_scale
, ic
->ewaldcoeff_q
, v_q_ewald_lr
);
1877 if (EVDW_PME(ic
->vdwtype
))
1879 snew_aligned(ic
->tabq_vdw_FDV0
, ic
->tabq_size
*4, 32);
1880 snew_aligned(ic
->tabq_vdw_F
, ic
->tabq_size
, 32);
1881 snew_aligned(ic
->tabq_vdw_V
, ic
->tabq_size
, 32);
1882 table_spline3_fill_ewald_lr(ic
->tabq_vdw_F
, ic
->tabq_vdw_V
, ic
->tabq_vdw_FDV0
,
1883 ic
->tabq_size
, 1/ic
->tabq_scale
, ic
->ewaldcoeff_lj
, v_lj_ewald_lr
);
1887 void init_interaction_const_tables(FILE *fp
,
1888 interaction_const_t
*ic
,
1891 if (EEL_PME_EWALD(ic
->eeltype
) || EVDW_PME(ic
->vdwtype
))
1893 init_ewald_f_table(ic
, rtab
);
1897 fprintf(fp
, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1898 1/ic
->tabq_scale
, ic
->tabq_size
);
1903 static void clear_force_switch_constants(shift_consts_t
*sc
)
1910 static void force_switch_constants(real p
,
1914 /* Here we determine the coefficient for shifting the force to zero
1915 * between distance rsw and the cut-off rc.
1916 * For a potential of r^-p, we have force p*r^-(p+1).
1917 * But to save flops we absorb p in the coefficient.
1919 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1920 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1922 sc
->c2
= ((p
+ 1)*rsw
- (p
+ 4)*rc
)/(pow(rc
, p
+ 2)*gmx::square(rc
- rsw
));
1923 sc
->c3
= -((p
+ 1)*rsw
- (p
+ 3)*rc
)/(pow(rc
, p
+ 2)*gmx::power3(rc
- rsw
));
1924 sc
->cpot
= -pow(rc
, -p
) + p
*sc
->c2
/3*gmx::power3(rc
- rsw
) + p
*sc
->c3
/4*gmx::power4(rc
- rsw
);
1927 static void potential_switch_constants(real rsw
, real rc
,
1928 switch_consts_t
*sc
)
1930 /* The switch function is 1 at rsw and 0 at rc.
1931 * The derivative and second derivate are zero at both ends.
1932 * rsw = max(r - r_switch, 0)
1933 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1934 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1935 * force = force*dsw - potential*sw
1938 sc
->c3
= -10/gmx::power3(rc
- rsw
);
1939 sc
->c4
= 15/gmx::power4(rc
- rsw
);
1940 sc
->c5
= -6/gmx::power5(rc
- rsw
);
1943 /*! \brief Construct interaction constants
1945 * This data is used (particularly) by search and force code for
1946 * short-range interactions. Many of these are constant for the whole
1947 * simulation; some are constant only after PME tuning completes.
1950 init_interaction_const(FILE *fp
,
1951 interaction_const_t
**interaction_const
,
1952 const t_inputrec
*ir
,
1953 const gmx_mtop_t
*mtop
,
1954 bool systemHasNetCharge
)
1956 interaction_const_t
*ic
;
1960 ic
->cutoff_scheme
= ir
->cutoff_scheme
;
1962 /* Just allocate something so we can free it */
1963 snew_aligned(ic
->tabq_coul_FDV0
, 16, 32);
1964 snew_aligned(ic
->tabq_coul_F
, 16, 32);
1965 snew_aligned(ic
->tabq_coul_V
, 16, 32);
1968 ic
->vdwtype
= ir
->vdwtype
;
1969 ic
->vdw_modifier
= ir
->vdw_modifier
;
1970 ic
->reppow
= mtop
->ffparams
.reppow
;
1971 ic
->rvdw
= cutoff_inf(ir
->rvdw
);
1972 ic
->rvdw_switch
= ir
->rvdw_switch
;
1973 ic
->ljpme_comb_rule
= ir
->ljpme_combination_rule
;
1974 ic
->useBuckingham
= (mtop
->ffparams
.functype
[0] == F_BHAM
);
1975 if (ic
->useBuckingham
)
1977 ic
->buckinghamBMax
= calcBuckinghamBMax(fp
, mtop
);
1980 initVdwEwaldParameters(fp
, ir
, ic
);
1982 clear_force_switch_constants(&ic
->dispersion_shift
);
1983 clear_force_switch_constants(&ic
->repulsion_shift
);
1985 switch (ic
->vdw_modifier
)
1987 case eintmodPOTSHIFT
:
1988 /* Only shift the potential, don't touch the force */
1989 ic
->dispersion_shift
.cpot
= -1.0/gmx::power6(ic
->rvdw
);
1990 ic
->repulsion_shift
.cpot
= -1.0/gmx::power12(ic
->rvdw
);
1992 case eintmodFORCESWITCH
:
1993 /* Switch the force, switch and shift the potential */
1994 force_switch_constants(6.0, ic
->rvdw_switch
, ic
->rvdw
,
1995 &ic
->dispersion_shift
);
1996 force_switch_constants(12.0, ic
->rvdw_switch
, ic
->rvdw
,
1997 &ic
->repulsion_shift
);
1999 case eintmodPOTSWITCH
:
2000 /* Switch the potential and force */
2001 potential_switch_constants(ic
->rvdw_switch
, ic
->rvdw
,
2005 case eintmodEXACTCUTOFF
:
2006 /* Nothing to do here */
2009 gmx_incons("unimplemented potential modifier");
2012 ic
->sh_invrc6
= -ic
->dispersion_shift
.cpot
;
2014 /* Electrostatics */
2015 ic
->eeltype
= ir
->coulombtype
;
2016 ic
->coulomb_modifier
= ir
->coulomb_modifier
;
2017 ic
->rcoulomb
= cutoff_inf(ir
->rcoulomb
);
2018 ic
->rcoulomb_switch
= ir
->rcoulomb_switch
;
2019 ic
->epsilon_r
= ir
->epsilon_r
;
2021 /* Set the Coulomb energy conversion factor */
2022 if (ic
->epsilon_r
!= 0)
2024 ic
->epsfac
= ONE_4PI_EPS0
/ic
->epsilon_r
;
2028 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2032 /* Reaction-field */
2033 if (EEL_RF(ic
->eeltype
))
2035 ic
->epsilon_rf
= ir
->epsilon_rf
;
2036 /* Generalized reaction field parameters are updated every step */
2037 if (ic
->eeltype
!= eelGRF
)
2039 calc_rffac(fp
, ic
->eeltype
, ic
->epsilon_r
, ic
->epsilon_rf
,
2040 ic
->rcoulomb
, 0, 0, nullptr,
2041 &ic
->k_rf
, &ic
->c_rf
);
2044 if (ir
->cutoff_scheme
== ecutsGROUP
&& ic
->eeltype
== eelRF_ZERO
)
2046 /* grompp should have done this, but this scheme is obsolete */
2047 ic
->coulomb_modifier
= eintmodEXACTCUTOFF
;
2052 /* For plain cut-off we might use the reaction-field kernels */
2053 ic
->epsilon_rf
= ic
->epsilon_r
;
2055 if (ir
->coulomb_modifier
== eintmodPOTSHIFT
)
2057 ic
->c_rf
= 1/ic
->rcoulomb
;
2065 initCoulombEwaldParameters(fp
, ir
, systemHasNetCharge
, ic
);
2069 real dispersion_shift
;
2071 dispersion_shift
= ic
->dispersion_shift
.cpot
;
2072 if (EVDW_PME(ic
->vdwtype
))
2074 dispersion_shift
-= ic
->sh_lj_ewald
;
2076 fprintf(fp
, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2077 ic
->repulsion_shift
.cpot
, dispersion_shift
);
2079 if (ic
->eeltype
== eelCUT
)
2081 fprintf(fp
, ", Coulomb %.e", -ic
->c_rf
);
2083 else if (EEL_PME(ic
->eeltype
))
2085 fprintf(fp
, ", Ewald %.3e", -ic
->sh_ewald
);
2090 *interaction_const
= ic
;
2094 done_interaction_const(interaction_const_t
*interaction_const
)
2096 sfree_aligned(interaction_const
->tabq_coul_FDV0
);
2097 sfree_aligned(interaction_const
->tabq_coul_F
);
2098 sfree_aligned(interaction_const
->tabq_coul_V
);
2099 sfree(interaction_const
);
2102 static void init_nb_verlet(const gmx::MDLogger
&mdlog
,
2103 nonbonded_verlet_t
**nb_verlet
,
2104 gmx_bool bFEP_NonBonded
,
2105 const t_inputrec
*ir
,
2106 const t_forcerec
*fr
,
2107 const t_commrec
*cr
,
2108 const gmx_hw_info_t
&hardwareInfo
,
2109 const gmx_device_info_t
*deviceInfo
,
2110 const gmx_mtop_t
*mtop
,
2113 nonbonded_verlet_t
*nbv
;
2116 nbv
= new nonbonded_verlet_t();
2118 nbv
->emulateGpu
= ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes
: EmulateGpuNonbonded::No
);
2119 nbv
->bUseGPU
= deviceInfo
!= nullptr;
2121 GMX_RELEASE_ASSERT(!(nbv
->emulateGpu
== EmulateGpuNonbonded::Yes
&& nbv
->bUseGPU
), "When GPU emulation is active, there cannot be a GPU assignment");
2124 nbv
->min_ci_balanced
= 0;
2126 nbv
->ngrp
= (DOMAINDECOMP(cr
) ? 2 : 1);
2127 for (int i
= 0; i
< nbv
->ngrp
; i
++)
2129 nbv
->grp
[i
].nbl_lists
.nnbl
= 0;
2130 nbv
->grp
[i
].kernel_type
= nbnxnkNotSet
;
2132 if (i
== 0) /* local */
2134 pick_nbnxn_kernel(mdlog
, fr
->use_simd_kernels
, hardwareInfo
,
2135 nbv
->bUseGPU
, nbv
->emulateGpu
, ir
,
2136 &nbv
->grp
[i
].kernel_type
,
2137 &nbv
->grp
[i
].ewald_excl
,
2140 else /* non-local */
2142 /* Use the same kernel for local and non-local interactions */
2143 nbv
->grp
[i
].kernel_type
= nbv
->grp
[0].kernel_type
;
2144 nbv
->grp
[i
].ewald_excl
= nbv
->grp
[0].ewald_excl
;
2148 nbv
->listParams
= gmx::compat::make_unique
<NbnxnListParameters
>(ir
->rlist
);
2149 setupDynamicPairlistPruning(mdlog
, ir
, mtop
, box
, nbv
->grp
[0].kernel_type
, fr
->ic
,
2150 nbv
->listParams
.get());
2152 nbv
->nbs
= gmx::compat::make_unique
<nbnxn_search
>(DOMAINDECOMP(cr
) ? &cr
->dd
->nc
: nullptr,
2153 DOMAINDECOMP(cr
) ? domdec_zones(cr
->dd
) : nullptr,
2155 gmx_omp_nthreads_get(emntPairsearch
));
2157 for (int i
= 0; i
< nbv
->ngrp
; i
++)
2159 nbnxn_init_pairlist_set(&nbv
->grp
[i
].nbl_lists
,
2160 nbnxn_kernel_pairlist_simple(nbv
->grp
[i
].kernel_type
),
2161 /* 8x8x8 "non-simple" lists are ATM always combined */
2162 !nbnxn_kernel_pairlist_simple(nbv
->grp
[i
].kernel_type
));
2165 int enbnxninitcombrule
;
2166 if (fr
->ic
->vdwtype
== evdwCUT
&&
2167 (fr
->ic
->vdw_modifier
== eintmodNONE
||
2168 fr
->ic
->vdw_modifier
== eintmodPOTSHIFT
) &&
2169 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2171 /* Plain LJ cut-off: we can optimize with combination rules */
2172 enbnxninitcombrule
= enbnxninitcombruleDETECT
;
2174 else if (fr
->ic
->vdwtype
== evdwPME
)
2176 /* LJ-PME: we need to use a combination rule for the grid */
2177 if (fr
->ljpme_combination_rule
== eljpmeGEOM
)
2179 enbnxninitcombrule
= enbnxninitcombruleGEOM
;
2183 enbnxninitcombrule
= enbnxninitcombruleLB
;
2188 /* We use a full combination matrix: no rule required */
2189 enbnxninitcombrule
= enbnxninitcombruleNONE
;
2192 nbv
->nbat
= new nbnxn_atomdata_t(nbv
->bUseGPU
? gmx::PinningPolicy::PinnedIfSupported
: gmx::PinningPolicy::CannotBePinned
);
2193 int mimimumNumEnergyGroupNonbonded
= ir
->opts
.ngener
;
2194 if (ir
->opts
.ngener
- ir
->nwall
== 1)
2196 /* We have only one non-wall energy group, we do not need energy group
2197 * support in the non-bondeds kernels, since all non-bonded energy
2198 * contributions go to the first element of the energy group matrix.
2200 mimimumNumEnergyGroupNonbonded
= 1;
2202 bool bSimpleList
= nbnxn_kernel_pairlist_simple(nbv
->grp
[0].kernel_type
);
2203 nbnxn_atomdata_init(mdlog
,
2205 nbv
->grp
[0].kernel_type
,
2207 fr
->ntype
, fr
->nbfp
,
2208 mimimumNumEnergyGroupNonbonded
,
2209 bSimpleList
? gmx_omp_nthreads_get(emntNonbonded
) : 1);
2213 /* init the NxN GPU data; the last argument tells whether we'll have
2214 * both local and non-local NB calculation on GPU */
2215 nbnxn_gpu_init(&nbv
->gpu_nbv
,
2218 nbv
->listParams
.get(),
2223 if ((env
= getenv("GMX_NB_MIN_CI")) != nullptr)
2227 nbv
->min_ci_balanced
= strtol(env
, &end
, 10);
2228 if (!end
|| (*end
!= 0) || nbv
->min_ci_balanced
< 0)
2230 gmx_fatal(FARGS
, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env
);
2235 fprintf(debug
, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2236 nbv
->min_ci_balanced
);
2241 nbv
->min_ci_balanced
= nbnxn_gpu_min_ci_balanced(nbv
->gpu_nbv
);
2244 fprintf(debug
, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2245 nbv
->min_ci_balanced
);
2254 gmx_bool
usingGpu(nonbonded_verlet_t
*nbv
)
2256 return nbv
!= nullptr && nbv
->bUseGPU
;
2259 void init_forcerec(FILE *fp
,
2260 const gmx::MDLogger
&mdlog
,
2263 const t_inputrec
*ir
,
2264 const gmx_mtop_t
*mtop
,
2265 const t_commrec
*cr
,
2269 gmx::ArrayRef
<const std::string
> tabbfnm
,
2270 const gmx_hw_info_t
&hardwareInfo
,
2271 const gmx_device_info_t
*deviceInfo
,
2272 const bool useGpuForBonded
,
2273 gmx_bool bNoSolvOpt
,
2276 int m
, negp_pp
, negptable
, egi
, egj
;
2281 gmx_bool bGenericKernelOnly
;
2282 gmx_bool needGroupSchemeTables
, bSomeNormalNbListsAreInUse
;
2283 gmx_bool bFEP_NonBonded
;
2284 int *nm_ind
, egp_flags
;
2286 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2287 fr
->use_simd_kernels
= TRUE
;
2289 if (check_box(ir
->ePBC
, box
))
2291 gmx_fatal(FARGS
, "%s", check_box(ir
->ePBC
, box
));
2294 /* Test particle insertion ? */
2297 /* Set to the size of the molecule to be inserted (the last one) */
2298 /* Because of old style topologies, we have to use the last cg
2299 * instead of the last molecule type.
2301 cgs
= &mtop
->moltype
[mtop
->molblock
.back().type
].cgs
;
2302 fr
->n_tpi
= cgs
->index
[cgs
->nr
] - cgs
->index
[cgs
->nr
-1];
2303 gmx::RangePartitioning molecules
= gmx_mtop_molecules(*mtop
);
2304 if (fr
->n_tpi
!= molecules
.block(molecules
.numBlocks() - 1).size())
2306 gmx_fatal(FARGS
, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2314 if (ir
->coulombtype
== eelRF_NEC_UNSUPPORTED
)
2316 gmx_fatal(FARGS
, "%s electrostatics is no longer supported",
2317 eel_names
[ir
->coulombtype
]);
2322 gmx_fatal(FARGS
, "AdResS simulations are no longer supported");
2324 if (ir
->useTwinRange
)
2326 gmx_fatal(FARGS
, "Twin-range simulations are no longer supported");
2328 /* Copy the user determined parameters */
2329 fr
->userint1
= ir
->userint1
;
2330 fr
->userint2
= ir
->userint2
;
2331 fr
->userint3
= ir
->userint3
;
2332 fr
->userint4
= ir
->userint4
;
2333 fr
->userreal1
= ir
->userreal1
;
2334 fr
->userreal2
= ir
->userreal2
;
2335 fr
->userreal3
= ir
->userreal3
;
2336 fr
->userreal4
= ir
->userreal4
;
2339 fr
->fc_stepsize
= ir
->fc_stepsize
;
2342 fr
->efep
= ir
->efep
;
2343 fr
->sc_alphavdw
= ir
->fepvals
->sc_alpha
;
2344 if (ir
->fepvals
->bScCoul
)
2346 fr
->sc_alphacoul
= ir
->fepvals
->sc_alpha
;
2347 fr
->sc_sigma6_min
= gmx::power6(ir
->fepvals
->sc_sigma_min
);
2351 fr
->sc_alphacoul
= 0;
2352 fr
->sc_sigma6_min
= 0; /* only needed when bScCoul is on */
2354 fr
->sc_power
= ir
->fepvals
->sc_power
;
2355 fr
->sc_r_power
= ir
->fepvals
->sc_r_power
;
2356 fr
->sc_sigma6_def
= gmx::power6(ir
->fepvals
->sc_sigma
);
2358 env
= getenv("GMX_SCSIGMA_MIN");
2362 sscanf(env
, "%20lf", &dbl
);
2363 fr
->sc_sigma6_min
= gmx::power6(dbl
);
2366 fprintf(fp
, "Setting the minimum soft core sigma to %g nm\n", dbl
);
2370 fr
->bNonbonded
= TRUE
;
2371 if (getenv("GMX_NO_NONBONDED") != nullptr)
2373 /* turn off non-bonded calculations */
2374 fr
->bNonbonded
= FALSE
;
2375 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
2376 "Found environment variable GMX_NO_NONBONDED.\n"
2377 "Disabling nonbonded calculations.");
2380 bGenericKernelOnly
= FALSE
;
2382 /* We now check in the NS code whether a particular combination of interactions
2383 * can be used with water optimization, and disable it if that is not the case.
2386 if (getenv("GMX_NB_GENERIC") != nullptr)
2391 "Found environment variable GMX_NB_GENERIC.\n"
2392 "Disabling all interaction-specific nonbonded kernels, will only\n"
2393 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2395 bGenericKernelOnly
= TRUE
;
2398 if (bGenericKernelOnly
)
2403 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
2405 fr
->use_simd_kernels
= FALSE
;
2409 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2410 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2411 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2415 fr
->bBHAM
= (mtop
->ffparams
.functype
[0] == F_BHAM
);
2417 /* Neighbour searching stuff */
2418 fr
->cutoff_scheme
= ir
->cutoff_scheme
;
2419 fr
->bGrid
= (ir
->ns_type
== ensGRID
);
2420 fr
->ePBC
= ir
->ePBC
;
2422 if (fr
->cutoff_scheme
== ecutsGROUP
)
2424 const char *note
= "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2425 "removed in a future release when 'verlet' supports all interaction forms.\n";
2429 fprintf(stderr
, "\n%s\n", note
);
2433 fprintf(fp
, "\n%s\n", note
);
2437 /* Determine if we will do PBC for distances in bonded interactions */
2438 if (fr
->ePBC
== epbcNONE
)
2440 fr
->bMolPBC
= FALSE
;
2444 if (!DOMAINDECOMP(cr
))
2448 bSHAKE
= (ir
->eConstrAlg
== econtSHAKE
&&
2449 (gmx_mtop_ftype_count(mtop
, F_CONSTR
) > 0 ||
2450 gmx_mtop_ftype_count(mtop
, F_CONSTRNC
) > 0));
2452 /* The group cut-off scheme and SHAKE assume charge groups
2453 * are whole, but not using molpbc is faster in most cases.
2454 * With intermolecular interactions we need PBC for calculating
2455 * distances between atoms in different molecules.
2457 if ((fr
->cutoff_scheme
== ecutsGROUP
|| bSHAKE
) &&
2458 !mtop
->bIntermolecularInteractions
)
2460 fr
->bMolPBC
= ir
->bPeriodicMols
;
2462 if (bSHAKE
&& fr
->bMolPBC
)
2464 gmx_fatal(FARGS
, "SHAKE is not supported with periodic molecules");
2469 /* Not making molecules whole is faster in most cases,
2470 * but With orientation restraints we need whole molecules.
2472 fr
->bMolPBC
= (fcd
->orires
.nr
== 0);
2474 if (getenv("GMX_USE_GRAPH") != nullptr)
2476 fr
->bMolPBC
= FALSE
;
2479 GMX_LOG(mdlog
.warning
).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
2482 if (mtop
->bIntermolecularInteractions
)
2484 GMX_LOG(mdlog
.warning
).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
2488 GMX_RELEASE_ASSERT(fr
->bMolPBC
|| !mtop
->bIntermolecularInteractions
, "We need to use PBC within molecules with inter-molecular interactions");
2490 if (bSHAKE
&& fr
->bMolPBC
)
2492 gmx_fatal(FARGS
, "SHAKE is not properly supported with intermolecular interactions. For short simulations where linked molecules remain in the same periodic image, the environment variable GMX_USE_GRAPH can be used to override this check.\n");
2498 fr
->bMolPBC
= dd_bonded_molpbc(cr
->dd
, fr
->ePBC
);
2502 fr
->rc_scaling
= ir
->refcoord_scaling
;
2503 copy_rvec(ir
->posres_com
, fr
->posres_com
);
2504 copy_rvec(ir
->posres_comB
, fr
->posres_comB
);
2505 fr
->rlist
= cutoff_inf(ir
->rlist
);
2506 fr
->ljpme_combination_rule
= ir
->ljpme_combination_rule
;
2508 /* This now calculates sum for q and c6*/
2509 bool systemHasNetCharge
= set_chargesum(fp
, fr
, mtop
);
2511 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2512 init_interaction_const(fp
, &fr
->ic
, ir
, mtop
, systemHasNetCharge
);
2513 init_interaction_const_tables(fp
, fr
->ic
, ir
->rlist
+ ir
->tabext
);
2515 const interaction_const_t
*ic
= fr
->ic
;
2517 /* TODO: Replace this Ewald table or move it into interaction_const_t */
2518 if (ir
->coulombtype
== eelEWALD
)
2520 init_ewald_tab(&(fr
->ewald_table
), ir
, fp
);
2523 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2524 switch (ic
->eeltype
)
2527 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_COULOMB
;
2532 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_REACTIONFIELD
;
2536 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_REACTIONFIELD
;
2537 GMX_RELEASE_ASSERT(ic
->coulomb_modifier
== eintmodEXACTCUTOFF
, "With the group scheme RF-zero needs the exact cut-off modifier");
2546 case eelPMEUSERSWITCH
:
2547 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_CUBICSPLINETABLE
;
2553 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_EWALD
;
2557 gmx_fatal(FARGS
, "Unsupported electrostatic interaction: %s", eel_names
[ic
->eeltype
]);
2559 fr
->nbkernel_elec_modifier
= ic
->coulomb_modifier
;
2561 /* Vdw: Translate from mdp settings to kernel format */
2562 switch (ic
->vdwtype
)
2567 fr
->nbkernel_vdw_interaction
= GMX_NBKERNEL_VDW_BUCKINGHAM
;
2571 fr
->nbkernel_vdw_interaction
= GMX_NBKERNEL_VDW_LENNARDJONES
;
2575 fr
->nbkernel_vdw_interaction
= GMX_NBKERNEL_VDW_LJEWALD
;
2581 case evdwENCADSHIFT
:
2582 fr
->nbkernel_vdw_interaction
= GMX_NBKERNEL_VDW_CUBICSPLINETABLE
;
2586 gmx_fatal(FARGS
, "Unsupported vdw interaction: %s", evdw_names
[ic
->vdwtype
]);
2588 fr
->nbkernel_vdw_modifier
= ic
->vdw_modifier
;
2590 if (ir
->cutoff_scheme
== ecutsGROUP
)
2592 fr
->bvdwtab
= ((ic
->vdwtype
!= evdwCUT
|| !gmx_within_tol(ic
->reppow
, 12.0, 10*GMX_DOUBLE_EPS
))
2593 && !EVDW_PME(ic
->vdwtype
));
2594 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2595 fr
->bcoultab
= !(ic
->eeltype
== eelCUT
||
2596 ic
->eeltype
== eelEWALD
||
2597 ic
->eeltype
== eelPME
||
2598 ic
->eeltype
== eelP3M_AD
||
2599 ic
->eeltype
== eelRF
||
2600 ic
->eeltype
== eelRF_ZERO
);
2602 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2603 * going to be faster to tabulate the interaction than calling the generic kernel.
2604 * However, if generic kernels have been requested we keep things analytically.
2606 if (fr
->nbkernel_elec_modifier
== eintmodPOTSWITCH
&&
2607 fr
->nbkernel_vdw_modifier
== eintmodPOTSWITCH
&&
2608 !bGenericKernelOnly
)
2610 if ((ic
->rcoulomb_switch
!= ic
->rvdw_switch
) || (ic
->rcoulomb
!= ic
->rvdw
))
2612 fr
->bcoultab
= TRUE
;
2613 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2614 * which would otherwise need two tables.
2618 else if ((fr
->nbkernel_elec_modifier
== eintmodPOTSHIFT
&& fr
->nbkernel_vdw_modifier
== eintmodPOTSHIFT
) ||
2619 ((fr
->nbkernel_elec_interaction
== GMX_NBKERNEL_ELEC_REACTIONFIELD
&&
2620 fr
->nbkernel_elec_modifier
== eintmodEXACTCUTOFF
&&
2621 (fr
->nbkernel_vdw_modifier
== eintmodPOTSWITCH
|| fr
->nbkernel_vdw_modifier
== eintmodPOTSHIFT
))))
2623 if ((ic
->rcoulomb
!= ic
->rvdw
) && (!bGenericKernelOnly
))
2625 fr
->bcoultab
= TRUE
;
2629 if (fr
->nbkernel_elec_modifier
== eintmodFORCESWITCH
)
2631 fr
->bcoultab
= TRUE
;
2633 if (fr
->nbkernel_vdw_modifier
== eintmodFORCESWITCH
)
2638 if (getenv("GMX_REQUIRE_TABLES"))
2641 fr
->bcoultab
= TRUE
;
2646 fprintf(fp
, "Table routines are used for coulomb: %s\n",
2647 gmx::boolToString(fr
->bcoultab
));
2648 fprintf(fp
, "Table routines are used for vdw: %s\n",
2649 gmx::boolToString(fr
->bvdwtab
));
2654 fr
->nbkernel_vdw_interaction
= GMX_NBKERNEL_VDW_CUBICSPLINETABLE
;
2655 fr
->nbkernel_vdw_modifier
= eintmodNONE
;
2659 fr
->nbkernel_elec_interaction
= GMX_NBKERNEL_ELEC_CUBICSPLINETABLE
;
2660 fr
->nbkernel_elec_modifier
= eintmodNONE
;
2664 if (ir
->cutoff_scheme
== ecutsVERLET
)
2666 if (!gmx_within_tol(ic
->reppow
, 12.0, 10*GMX_DOUBLE_EPS
))
2668 gmx_fatal(FARGS
, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names
[ir
->cutoff_scheme
]);
2670 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2671 * while mdrun does not (and never did) support this.
2673 if (EEL_USER(fr
->ic
->eeltype
))
2675 gmx_fatal(FARGS
, "Combination of %s and cutoff scheme %s is not supported",
2676 eel_names
[ir
->coulombtype
], ecutscheme_names
[ir
->cutoff_scheme
]);
2679 fr
->bvdwtab
= FALSE
;
2680 fr
->bcoultab
= FALSE
;
2683 /* 1-4 interaction electrostatics */
2684 fr
->fudgeQQ
= mtop
->ffparams
.fudgeQQ
;
2686 /* Parameters for generalized RF */
2690 if (ic
->eeltype
== eelGRF
)
2692 init_generalized_rf(fp
, mtop
, ir
, fr
);
2695 fr
->haveDirectVirialContributions
=
2696 (EEL_FULL(ic
->eeltype
) || EVDW_PME(ic
->vdwtype
) ||
2697 fr
->forceProviders
->hasForceProvider() ||
2698 gmx_mtop_ftype_count(mtop
, F_POSRES
) > 0 ||
2699 gmx_mtop_ftype_count(mtop
, F_FBPOSRES
) > 0 ||
2705 if (fr
->haveDirectVirialContributions
)
2707 fr
->forceBufferForDirectVirialContributions
= new std::vector
<gmx::RVec
>;
2710 if (fr
->cutoff_scheme
== ecutsGROUP
&&
2711 ncg_mtop(mtop
) > fr
->cg_nalloc
&& !DOMAINDECOMP(cr
))
2713 /* Count the total number of charge groups */
2714 fr
->cg_nalloc
= ncg_mtop(mtop
);
2715 srenew(fr
->cg_cm
, fr
->cg_nalloc
);
2717 if (fr
->shift_vec
== nullptr)
2719 snew(fr
->shift_vec
, SHIFTS
);
2722 if (fr
->fshift
== nullptr)
2724 snew(fr
->fshift
, SHIFTS
);
2727 if (fr
->nbfp
== nullptr)
2729 fr
->ntype
= mtop
->ffparams
.atnr
;
2730 fr
->nbfp
= mk_nbfp(&mtop
->ffparams
, fr
->bBHAM
);
2731 if (EVDW_PME(ic
->vdwtype
))
2733 fr
->ljpme_c6grid
= make_ljpme_c6grid(&mtop
->ffparams
, fr
);
2737 /* Copy the energy group exclusions */
2738 fr
->egp_flags
= ir
->opts
.egp_flags
;
2740 /* Van der Waals stuff */
2741 if ((ic
->vdwtype
!= evdwCUT
) && (ic
->vdwtype
!= evdwUSER
) && !fr
->bBHAM
)
2743 if (ic
->rvdw_switch
>= ic
->rvdw
)
2745 gmx_fatal(FARGS
, "rvdw_switch (%f) must be < rvdw (%f)",
2746 ic
->rvdw_switch
, ic
->rvdw
);
2750 fprintf(fp
, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2751 (ic
->eeltype
== eelSWITCH
) ? "switched" : "shifted",
2752 ic
->rvdw_switch
, ic
->rvdw
);
2756 if (fr
->bBHAM
&& EVDW_PME(ic
->vdwtype
))
2758 gmx_fatal(FARGS
, "LJ PME not supported with Buckingham");
2761 if (fr
->bBHAM
&& (ic
->vdwtype
== evdwSHIFT
|| ic
->vdwtype
== evdwSWITCH
))
2763 gmx_fatal(FARGS
, "Switch/shift interaction not supported with Buckingham");
2766 if (fr
->bBHAM
&& fr
->cutoff_scheme
== ecutsVERLET
)
2768 gmx_fatal(FARGS
, "Verlet cutoff-scheme is not supported with Buckingham");
2771 if (fp
&& fr
->cutoff_scheme
== ecutsGROUP
)
2773 fprintf(fp
, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2774 fr
->rlist
, ic
->rcoulomb
, fr
->bBHAM
? "BHAM" : "LJ", ic
->rvdw
);
2777 fr
->eDispCorr
= ir
->eDispCorr
;
2778 fr
->numAtomsForDispersionCorrection
= mtop
->natoms
;
2779 if (ir
->eDispCorr
!= edispcNO
)
2781 set_avcsixtwelve(fp
, fr
, mtop
);
2784 if (ir
->implicit_solvent
)
2786 gmx_fatal(FARGS
, "Implict solvation is no longer supported.");
2789 /* Construct tables for the group scheme. A little unnecessary to
2790 * make both vdw and coul tables sometimes, but what the
2791 * heck. Note that both cutoff schemes construct Ewald tables in
2792 * init_interaction_const_tables. */
2793 needGroupSchemeTables
= (ir
->cutoff_scheme
== ecutsGROUP
&&
2794 (fr
->bcoultab
|| fr
->bvdwtab
));
2796 negp_pp
= ir
->opts
.ngener
- ir
->nwall
;
2798 if (!needGroupSchemeTables
)
2800 bSomeNormalNbListsAreInUse
= TRUE
;
2805 bSomeNormalNbListsAreInUse
= FALSE
;
2806 for (egi
= 0; egi
< negp_pp
; egi
++)
2808 for (egj
= egi
; egj
< negp_pp
; egj
++)
2810 egp_flags
= ir
->opts
.egp_flags
[GID(egi
, egj
, ir
->opts
.ngener
)];
2811 if (!(egp_flags
& EGP_EXCL
))
2813 if (egp_flags
& EGP_TABLE
)
2819 bSomeNormalNbListsAreInUse
= TRUE
;
2824 if (bSomeNormalNbListsAreInUse
)
2826 fr
->nnblists
= negptable
+ 1;
2830 fr
->nnblists
= negptable
;
2832 if (fr
->nnblists
> 1)
2834 snew(fr
->gid2nblists
, ir
->opts
.ngener
*ir
->opts
.ngener
);
2838 snew(fr
->nblists
, fr
->nnblists
);
2840 /* This code automatically gives table length tabext without cut-off's,
2841 * in that case grompp should already have checked that we do not need
2842 * normal tables and we only generate tables for 1-4 interactions.
2844 rtab
= ir
->rlist
+ ir
->tabext
;
2846 if (needGroupSchemeTables
)
2848 /* make tables for ordinary interactions */
2849 if (bSomeNormalNbListsAreInUse
)
2851 make_nbf_tables(fp
, ic
, rtab
, tabfn
, nullptr, nullptr, &fr
->nblists
[0]);
2860 /* Read the special tables for certain energy group pairs */
2861 nm_ind
= mtop
->groups
.grps
[egcENER
].nm_ind
;
2862 for (egi
= 0; egi
< negp_pp
; egi
++)
2864 for (egj
= egi
; egj
< negp_pp
; egj
++)
2866 egp_flags
= ir
->opts
.egp_flags
[GID(egi
, egj
, ir
->opts
.ngener
)];
2867 if ((egp_flags
& EGP_TABLE
) && !(egp_flags
& EGP_EXCL
))
2869 if (fr
->nnblists
> 1)
2871 fr
->gid2nblists
[GID(egi
, egj
, ir
->opts
.ngener
)] = m
;
2873 /* Read the table file with the two energy groups names appended */
2874 make_nbf_tables(fp
, ic
, rtab
, tabfn
,
2875 *mtop
->groups
.grpname
[nm_ind
[egi
]],
2876 *mtop
->groups
.grpname
[nm_ind
[egj
]],
2880 else if (fr
->nnblists
> 1)
2882 fr
->gid2nblists
[GID(egi
, egj
, ir
->opts
.ngener
)] = 0;
2889 /* Tables might not be used for the potential modifier
2890 * interactions per se, but we still need them to evaluate
2891 * switch/shift dispersion corrections in this case. */
2892 if (fr
->eDispCorr
!= edispcNO
)
2894 fr
->dispersionCorrectionTable
= makeDispersionCorrectionTable(fp
, ic
, rtab
, tabfn
);
2897 /* We want to use unmodified tables for 1-4 coulombic
2898 * interactions, so we must in general have an extra set of
2900 if (gmx_mtop_ftype_count(mtop
, F_LJ14
) > 0 ||
2901 gmx_mtop_ftype_count(mtop
, F_LJC14_Q
) > 0 ||
2902 gmx_mtop_ftype_count(mtop
, F_LJC_PAIRS_NB
) > 0)
2904 fr
->pairsTable
= make_tables(fp
, ic
, tabpfn
, rtab
,
2905 GMX_MAKETABLES_14ONLY
);
2909 fr
->nwall
= ir
->nwall
;
2910 if (ir
->nwall
&& ir
->wall_type
== ewtTABLE
)
2912 make_wall_tables(fp
, ir
, tabfn
, &mtop
->groups
, fr
);
2915 if (fcd
&& !tabbfnm
.empty())
2917 // Need to catch std::bad_alloc
2918 // TODO Don't need to catch this here, when merging with master branch
2921 fcd
->bondtab
= make_bonded_tables(fp
,
2922 F_TABBONDS
, F_TABBONDSNC
,
2923 mtop
, tabbfnm
, "b");
2924 fcd
->angletab
= make_bonded_tables(fp
,
2926 mtop
, tabbfnm
, "a");
2927 fcd
->dihtab
= make_bonded_tables(fp
,
2929 mtop
, tabbfnm
, "d");
2931 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
2937 fprintf(debug
, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
2941 // QM/MM initialization if requested
2942 fr
->bQMMM
= ir
->bQMMM
;
2945 // Initialize QM/MM if supported
2948 GMX_LOG(mdlog
.info
).asParagraph().
2949 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
2950 "version. Please get in touch with the developers if you find the support useful, "
2951 "as help is needed if the functionality is to continue to be available.");
2952 fr
->qr
= mk_QMMMrec();
2953 init_QMMMrec(cr
, mtop
, ir
, fr
);
2957 gmx_incons("QM/MM was requested, but is only available when GROMACS "
2958 "is configured with QM/MM support");
2962 /* Set all the static charge group info */
2963 fr
->cginfo_mb
= init_cginfo_mb(fp
, mtop
, fr
, bNoSolvOpt
,
2965 &fr
->bExcl_IntraCGAll_InterCGNone
);
2966 if (DOMAINDECOMP(cr
))
2968 fr
->cginfo
= nullptr;
2972 fr
->cginfo
= cginfo_expand(mtop
->molblock
.size(), fr
->cginfo_mb
);
2975 if (!DOMAINDECOMP(cr
))
2977 forcerec_set_ranges(fr
, ncg_mtop(mtop
), ncg_mtop(mtop
),
2978 mtop
->natoms
, mtop
->natoms
, mtop
->natoms
);
2981 fr
->print_force
= print_force
;
2984 /* coarse load balancing vars */
2989 /* Initialize neighbor search */
2991 init_ns(fp
, cr
, fr
->ns
, fr
, mtop
);
2993 if (thisRankHasDuty(cr
, DUTY_PP
))
2995 gmx_nonbonded_setup(fr
, bGenericKernelOnly
);
2998 /* Initialize the thread working data for bonded interactions */
2999 init_bonded_threading(fp
, mtop
->groups
.grps
[egcENER
].nr
,
3000 &fr
->bondedThreading
);
3002 fr
->nthread_ewc
= gmx_omp_nthreads_get(emntBonded
);
3003 snew(fr
->ewc_t
, fr
->nthread_ewc
);
3005 if (fr
->cutoff_scheme
== ecutsVERLET
)
3007 // We checked the cut-offs in grompp, but double-check here.
3008 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3009 if (EEL_PME_EWALD(ir
->coulombtype
) && ir
->vdwtype
== eelCUT
)
3011 GMX_RELEASE_ASSERT(ir
->rcoulomb
>= ir
->rvdw
, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3015 GMX_RELEASE_ASSERT(ir
->rcoulomb
== ir
->rvdw
, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3018 init_nb_verlet(mdlog
, &fr
->nbv
, bFEP_NonBonded
, ir
, fr
,
3019 cr
, hardwareInfo
, deviceInfo
,
3022 if (useGpuForBonded
)
3024 auto stream
= DOMAINDECOMP(cr
) ?
3025 nbnxn_gpu_get_command_stream(fr
->nbv
->gpu_nbv
, eintNonlocal
) :
3026 nbnxn_gpu_get_command_stream(fr
->nbv
->gpu_nbv
, eintLocal
);
3027 // TODO the heap allocation is only needed while
3028 // t_forcerec lacks a constructor.
3029 fr
->gpuBonded
= new gmx::GpuBonded(mtop
->ffparams
,
3036 /* Here we switch from using mdlog, which prints the newline before
3037 * the paragraph, to our old fprintf logging, which prints the newline
3038 * after the paragraph, so we should add a newline here.
3043 if (ir
->eDispCorr
!= edispcNO
)
3045 calc_enervirdiff(fp
, ir
->eDispCorr
, fr
);
3049 /* Frees GPU memory and sets a tMPI node barrier.
3051 * Note that this function needs to be called even if GPUs are not used
3052 * in this run because the PME ranks have no knowledge of whether GPUs
3053 * are used or not, but all ranks need to enter the barrier below.
3054 * \todo Remove physical node barrier from this function after making sure
3055 * that it's not needed anymore (with a shared GPU run).
3057 void free_gpu_resources(t_forcerec
*fr
,
3058 const gmx::PhysicalNodeCommunicator
&physicalNodeCommunicator
)
3060 bool isPPrankUsingGPU
= (fr
!= nullptr) && (fr
->nbv
!= nullptr) && fr
->nbv
->bUseGPU
;
3062 /* stop the GPU profiler (only CUDA) */
3065 if (isPPrankUsingGPU
)
3067 /* free nbnxn data in GPU memory */
3068 nbnxn_gpu_free(fr
->nbv
->gpu_nbv
);
3069 delete fr
->gpuBonded
;
3070 fr
->gpuBonded
= nullptr;
3073 /* With tMPI we need to wait for all ranks to finish deallocation before
3074 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3077 * This is not a concern in OpenCL where we use one context per rank which
3078 * is freed in nbnxn_gpu_free().
3080 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
3081 * but it is easier and more futureproof to call it on the whole node.
3085 physicalNodeCommunicator
.barrier();
3089 void done_forcerec(t_forcerec
*fr
, int numMolBlocks
, int numEnergyGroups
)
3093 // PME-only ranks don't have a forcerec
3096 // cginfo is dynamically allocated if no domain decomposition
3097 if (fr
->cginfo
!= nullptr)
3101 done_cginfo_mb(fr
->cginfo_mb
, numMolBlocks
);
3103 done_interaction_const(fr
->ic
);
3104 sfree(fr
->shift_vec
);
3107 done_ns(fr
->ns
, numEnergyGroups
);
3109 tear_down_bonded_threading(fr
->bondedThreading
);
3110 GMX_RELEASE_ASSERT(fr
->gpuBonded
== nullptr, "Should have been deleted earlier, when used");
3111 fr
->bondedThreading
= nullptr;