Separate management of GPU contexts from modules
[gromacs.git] / src / gromacs / mdlib / forcerec.h
blob9906a755a118d0eaf068a11dcff4a947004dfed6
1 /*
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, 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/genborn.h"
42 #include "gromacs/mdlib/tgroup.h"
43 #include "gromacs/mdlib/vsite.h"
44 #include "gromacs/mdtypes/forcerec.h"
45 #include "gromacs/timing/wallcycle.h"
47 struct gmx_device_info_t;
48 struct t_commrec;
49 struct t_fcdata;
50 struct t_filenm;
52 namespace gmx
54 class MDLogger;
57 /*! \brief Create a new forcerec structure */
58 t_forcerec *mk_forcerec(void);
60 /*! \brief Print the contents of the forcerec to a file
62 * \param[in] fplog The log file to print to
63 * \param[in] fr The forcerec structure
65 void pr_forcerec(FILE *fplog, t_forcerec *fr);
67 /*! \brief Set the number of charge groups and atoms.
69 * The force calculation needs information on which atoms it
70 * should do work.
71 * \param[inout] fr The forcerec
72 * \param[in] ncg_home Number of charge groups on this processor
73 * \param[in] ncg_force Number of charge groups to compute force on
74 * \param[in] natoms_force Number of atoms to compute force on
75 * \param[in] natoms_force_constr Number of atoms involved in constraints
76 * \param[in] natoms_f_novirsum Number of atoms for which
77 * force is to be compute but no virial
79 void
80 forcerec_set_ranges(t_forcerec *fr,
81 int ncg_home, int ncg_force,
82 int natoms_force,
83 int natoms_force_constr, int natoms_f_novirsum);
85 /*! \brief Initiate table constants
87 * Initializes the tables in the interaction constant data structure.
88 * \param[in] fp File for debugging output
89 * \param[in] ic Structure holding the table constant
90 * \param[in] rtab The additional distance to add to tables
92 void init_interaction_const_tables(FILE *fp,
93 interaction_const_t *ic,
94 real rtab);
96 /*! \brief Initialize forcerec structure.
98 * The Force rec struct must be created with mk_forcerec.
99 * \param[in] fplog File for printing
100 * \param[in] mdlog File for printing
101 * \param[out] fr The forcerec
102 * \param[in] fcd Force constant data
103 * \param[in] ir Inputrec structure
104 * \param[in] mtop Molecular topology
105 * \param[in] cr Communication structures
106 * \param[in] box Simulation box
107 * \param[in] tabfn Table potential file for non-bonded interactions
108 * \param[in] tabpfn Table potential file for pair interactions
109 * \param[in] tabbfnm Table potential files for bonded interactions
110 * \param[in] deviceInfo Info about GPU device to use for short-ranged work
111 * \param[in] bNoSolvOpt Do not use solvent optimization
112 * \param[in] print_force Print forces for atoms with force >= print_force
114 void init_forcerec(FILE *fplog,
115 const gmx::MDLogger &mdlog,
116 t_forcerec *fr,
117 t_fcdata *fcd,
118 const t_inputrec *ir,
119 const gmx_mtop_t *mtop,
120 const t_commrec *cr,
121 matrix box,
122 const char *tabfn,
123 const char *tabpfn,
124 const t_filenm *tabbfnm,
125 const gmx_device_info_t *deviceInfo,
126 gmx_bool bNoSolvOpt,
127 real print_force);
129 /*! \brief Divide exclusions over threads
131 * Set the exclusion load for the local exclusions and possibly threads
132 * \param[out] fr The force record
133 * \param[in] top The topology
135 void forcerec_set_excl_load(t_forcerec *fr,
136 const gmx_localtop_t *top);
138 /*! \brief Update parameters dependent on box
140 * Updates parameters in the forcerec that are time dependent
141 * \param[out] fr The force record
142 * \param[in] box The simulation box
144 void update_forcerec(t_forcerec *fr, matrix box);
146 #endif