4d36ba7395ce2551d90d28e13c5b829e885f72ec
[alexxy/gromacs.git] / src / gromacs / mdlib / forcerec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifndef GMX_MDLIB_FORCEREC_H
38 #define GMX_MDLIB_FORCEREC_H
39
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"
46
47 struct gmx_device_info_t;
48 struct gmx_hw_info_t;
49 struct t_commrec;
50 struct t_fcdata;
51 struct t_filenm;
52 struct t_inputrec;
53
54 namespace gmx
55 {
56 class MDLogger;
57 class PhysicalNodeCommunicator;
58 }
59
60 /*! \brief Create a new forcerec structure */
61 t_forcerec *mk_forcerec(void);
62
63 //! Destroy a forcerec.
64 void done_forcerec(t_forcerec *fr, int numMolBlocks, int numEnergyGroups);
65
66 /*! \brief Print the contents of the forcerec to a file
67  *
68  * \param[in] fplog The log file to print to
69  * \param[in] fr    The forcerec structure
70  */
71 void pr_forcerec(FILE *fplog, t_forcerec *fr);
72
73 /*! \brief Set the number of charge groups and atoms.
74  *
75  * The force calculation needs information on which atoms it
76  * should do work.
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
84  */
85 void
86 forcerec_set_ranges(t_forcerec *fr,
87                     int ncg_home, int ncg_force,
88                     int natoms_force,
89                     int natoms_force_constr, int natoms_f_novirsum);
90
91 /*! \brief Initiate table constants
92  *
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
97  */
98 void init_interaction_const_tables(FILE                   *fp,
99                                    interaction_const_t    *ic,
100                                    real                    rtab);
101
102 /*! \brief Initialize forcerec structure.
103  *
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]  bNoSolvOpt  Do not use solvent optimization
119  * \param[in]  print_force Print forces for atoms with force >= print_force
120  */
121 void init_forcerec(FILE                             *fplog,
122                    const gmx::MDLogger              &mdlog,
123                    t_forcerec                       *fr,
124                    t_fcdata                         *fcd,
125                    const t_inputrec                 *ir,
126                    const gmx_mtop_t                 *mtop,
127                    const t_commrec                  *cr,
128                    matrix                            box,
129                    const char                       *tabfn,
130                    const char                       *tabpfn,
131                    gmx::ArrayRef<const std::string>  tabbfnm,
132                    const gmx_hw_info_t              &hardwareInfo,
133                    const gmx_device_info_t          *deviceInfo,
134                    gmx_bool                          bNoSolvOpt,
135                    real                              print_force);
136
137 /*! \brief Divide exclusions over threads
138  *
139  * Set the exclusion load for the local exclusions and possibly threads
140  * \param[out] fr  The force record
141  * \param[in]  top The topology
142  */
143 void forcerec_set_excl_load(t_forcerec           *fr,
144                             const gmx_localtop_t *top);
145
146 /*! \brief Update parameters dependent on box
147  *
148  * Updates parameters in the forcerec that are time dependent
149  * \param[out] fr  The force record
150  * \param[in]  box The simulation box
151  */
152 void update_forcerec(t_forcerec *fr, matrix box);
153
154 gmx_bool uses_simple_tables(int                 cutoff_scheme,
155                             nonbonded_verlet_t *nbv,
156                             int                 group);
157 /* Returns whether simple tables (i.e. not for use with GPUs) are used
158  * with the type of kernel indicated.
159  */
160
161 gmx_bool can_use_allvsall(const t_inputrec *ir,
162                           gmx_bool bPrintNote, const t_commrec *cr, FILE *fp);
163 /* Returns if we can use all-vs-all loops.
164  * If bPrintNote==TRUE, prints a note, if necessary, to stderr
165  * and fp (if !=NULL) on the master node.
166  */
167
168 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
169                               const t_inputrec    *ir);
170 /* Return if CPU SIMD support exists for the given inputrec
171  * If the return value is FALSE and fplog/cr != NULL, prints a fallback
172  * message to fplog/stderr.
173  */
174
175 /* Compute the average C6 and C12 params for LJ corrections */
176 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr,
177                       const gmx_mtop_t *mtop);
178
179 void free_gpu_resources(const t_forcerec                    *fr,
180                         const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator);
181
182 #endif