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.
37 #ifndef GMX_MDLIB_FORCEREC_H
38 #define GMX_MDLIB_FORCEREC_H
40 #include "gromacs/mdlib/force_flags.h"
41 #include "gromacs/mdlib/tgroup.h"
42 #include "gromacs/mdlib/vsite.h"
43 #include "gromacs/mdtypes/forcerec.h"
44 #include "gromacs/timing/wallcycle.h"
45 #include "gromacs/utility/arrayref.h"
47 struct gmx_device_info_t;
57 class PhysicalNodeCommunicator;
60 /*! \brief Create a new forcerec structure */
61 t_forcerec *mk_forcerec();
63 //! Destroy a forcerec.
64 void done_forcerec(t_forcerec *fr, int numMolBlocks, int numEnergyGroups);
66 /*! \brief Print the contents of the forcerec to a file
68 * \param[in] fplog The log file to print to
69 * \param[in] fr The forcerec structure
71 void pr_forcerec(FILE *fplog, t_forcerec *fr);
73 /*! \brief Set the number of charge groups and atoms.
75 * The force calculation needs information on which atoms it
77 * \param[inout] fr The forcerec
78 * \param[in] ncg_home Number of charge groups on this processor
79 * \param[in] ncg_force Number of charge groups to compute force on
80 * \param[in] natoms_force Number of atoms to compute force on
81 * \param[in] natoms_force_constr Number of atoms involved in constraints
82 * \param[in] natoms_f_novirsum Number of atoms for which
83 * force is to be compute but no virial
86 forcerec_set_ranges(t_forcerec *fr,
87 int ncg_home, int ncg_force,
89 int natoms_force_constr, int natoms_f_novirsum);
91 /*! \brief Initiate table constants
93 * Initializes the tables in the interaction constant data structure.
94 * \param[in] fp File for debugging output
95 * \param[in] ic Structure holding the table constant
96 * \param[in] rtab The additional distance to add to tables
98 void init_interaction_const_tables(FILE *fp,
99 interaction_const_t *ic,
102 /*! \brief Initialize forcerec structure.
104 * The Force rec struct must be created with mk_forcerec.
105 * \param[in] fplog File for printing
106 * \param[in] mdlog File for printing
107 * \param[out] fr The forcerec
108 * \param[in] fcd Force constant data
109 * \param[in] ir Inputrec structure
110 * \param[in] mtop Molecular topology
111 * \param[in] cr Communication structures
112 * \param[in] box Simulation box
113 * \param[in] tabfn Table potential file for non-bonded interactions
114 * \param[in] tabpfn Table potential file for pair interactions
115 * \param[in] tabbfnm Table potential files for bonded interactions
116 * \param[in] hardwareInfo Information about hardware
117 * \param[in] deviceInfo Info about GPU device to use for short-ranged work
118 * \param[in] useGpuForBonded Whether bonded interactions will run on a GPU
119 * \param[in] bNoSolvOpt Do not use solvent optimization
120 * \param[in] print_force Print forces for atoms with force >= print_force
122 void init_forcerec(FILE *fplog,
123 const gmx::MDLogger &mdlog,
126 const t_inputrec *ir,
127 const gmx_mtop_t *mtop,
132 gmx::ArrayRef<const std::string> tabbfnm,
133 const gmx_hw_info_t &hardwareInfo,
134 const gmx_device_info_t *deviceInfo,
135 bool useGpuForBonded,
139 /*! \brief Divide exclusions over threads
141 * Set the exclusion load for the local exclusions and possibly threads
142 * \param[out] fr The force record
143 * \param[in] top The topology
145 void forcerec_set_excl_load(t_forcerec *fr,
146 const gmx_localtop_t *top);
148 /*! \brief Update parameters dependent on box
150 * Updates parameters in the forcerec that are time dependent
151 * \param[out] fr The force record
152 * \param[in] box The simulation box
154 void update_forcerec(t_forcerec *fr, matrix box);
156 gmx_bool uses_simple_tables(int cutoff_scheme,
157 nonbonded_verlet_t *nbv,
159 /* Returns whether simple tables (i.e. not for use with GPUs) are used
160 * with the type of kernel indicated.
163 /* Compute the average C6 and C12 params for LJ corrections */
164 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr,
165 const gmx_mtop_t *mtop);
167 void free_gpu_resources(t_forcerec *fr,
168 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator);