Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / forcerec.h
blobaed645f1867def785f18a22fd570c5953069ca4b
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 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_MDLIB_FORCEREC_H
39 #define GMX_MDLIB_FORCEREC_H
41 #include "gromacs/math/vec.h"
42 #include "gromacs/timing/wallcycle.h"
43 #include "gromacs/utility/arrayref.h"
45 struct gmx_hw_info_t;
46 struct t_commrec;
47 struct t_forcerec;
48 struct t_filenm;
49 struct t_inputrec;
50 struct gmx_localtop_t;
51 struct gmx_mtop_t;
52 struct gmx_wallcycle;
53 struct interaction_const_t;
55 namespace gmx
57 class MDLogger;
58 class PhysicalNodeCommunicator;
59 } // namespace gmx
61 /*! \brief Print the contents of the forcerec to a file
63 * \param[in] fplog The log file to print to
64 * \param[in] fr The forcerec structure
66 void pr_forcerec(FILE* fplog, t_forcerec* fr);
68 /*! \brief Set the number of charge groups and atoms.
70 * The force calculation needs information on which atoms it
71 * should do work.
72 * \param[inout] fr The forcerec
73 * \param[in] natoms_force Number of atoms to compute force on
74 * \param[in] natoms_force_constr Number of atoms involved in constraints
75 * \param[in] natoms_f_novirsum Number of atoms for which
76 * force is to be compute but no virial
78 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum);
80 /*! \brief Initiate table constants
82 * Initializes the tables in the interaction constant data structure.
83 * \param[in] fp File for debugging output
84 * \param[in] ic Structure holding the table constant
85 * \param[in] tableExtensionLength Length by which to extend the tables. Taken from the input record.
87 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, real tableExtensionLength);
89 /*! \brief Initialize forcerec structure.
91 * \param[in] fplog File for printing
92 * \param[in] mdlog File for printing
93 * \param[out] fr The forcerec
94 * \param[in] ir Inputrec structure
95 * \param[in] mtop Molecular topology
96 * \param[in] cr Communication structures
97 * \param[in] box Simulation box
98 * \param[in] tabfn Table potential file for non-bonded interactions
99 * \param[in] tabpfn Table potential file for pair interactions
100 * \param[in] tabbfnm Table potential files for bonded interactions
101 * \param[in] print_force Print forces for atoms with force >= print_force
103 void init_forcerec(FILE* fplog,
104 const gmx::MDLogger& mdlog,
105 t_forcerec* fr,
106 const t_inputrec* ir,
107 const gmx_mtop_t* mtop,
108 const t_commrec* cr,
109 matrix box,
110 const char* tabfn,
111 const char* tabpfn,
112 gmx::ArrayRef<const std::string> tabbfnm,
113 real print_force);
115 /*! \brief Check whether molecules are ever distributed over PBC boundaries
117 * Note: This covers only the non-DD case. For DD runs, domdec.h offers an
118 * equivalent dd_bonded_molpbc(...) function.
120 * \param[in] ir Inputrec structure
121 * \param[in] mtop Molecular topology
122 * \param[in] mdlog File for printing
124 bool areMoleculesDistributedOverPbc(const t_inputrec& ir, const gmx_mtop_t& mtop, const gmx::MDLogger& mdlog);
126 /*! \brief Divide exclusions over threads
128 * Set the exclusion load for the local exclusions and possibly threads
129 * \param[out] fr The force record
130 * \param[in] top The topology
132 void forcerec_set_excl_load(t_forcerec* fr, const gmx_localtop_t* top);
134 #endif