Removed include of types/ifunc.h from typedefs.h
[gromacs.git] / src / gromacs / mdlib / forcerec.h
blob46d2eab05cfb06514c1597676f0bae1697d44528
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, 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/legacyheaders/genborn.h"
41 #include "gromacs/legacyheaders/network.h"
42 #include "gromacs/legacyheaders/tgroup.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/vsite.h"
45 #include "gromacs/legacyheaders/types/force_flags.h"
46 #include "gromacs/timing/wallcycle.h"
48 struct t_fcdata;
50 /*! \brief Create a new forcerec structure */
51 t_forcerec *mk_forcerec(void);
53 /*! \brief Print the contents of the forcerec to a file
55 * \param[in] fplog The log file to print to
56 * \param[in] fr The forcerec structure
58 void pr_forcerec(FILE *fplog, t_forcerec *fr);
60 /*! \brief Set the number of charge groups and atoms.
62 * The force calculation needs information on which atoms it
63 * should do work.
64 * \param[inout] fr The forcerec
65 * \param[in] ncg_home Number of charge groups on this processor
66 * \param[in] ncg_force Number of charge groups to compute force on
67 * \param[in] natoms_force Number of atoms to compute force on
68 * \param[in] natoms_force_constr Number of atoms involved in constraints
69 * \param[in] natoms_f_novirsum Number of atoms for which
70 * force is to be compute but no virial
72 void
73 forcerec_set_ranges(t_forcerec *fr,
74 int ncg_home, int ncg_force,
75 int natoms_force,
76 int natoms_force_constr, int natoms_f_novirsum);
78 /*! \brief Initiate table constants
80 * Initializes the tables in the interaction constant data structure.
81 * \param[in] fp File for debugging output
82 * \param[in] ic Structure holding the table constant
83 * \param[in] rtab The additional distance to add to tables
85 void init_interaction_const_tables(FILE *fp,
86 interaction_const_t *ic,
87 real rtab);
89 /*! \brief Initialize forcerec structure.
91 * The Force rec struct must be created with mk_forcerec.
92 * \param[in] fplog File for printing
93 * \param[in] oenv Output environment structure
94 * \param[out] fr The forcerec
95 * \param[in] fcd Force constant data
96 * \param[in] ir Inputrec structure
97 * \param[in] mtop Molecular topology
98 * \param[in] cr Communication structures
99 * \param[in] box Simulation box
100 * \param[in] tabfn Table potential file
101 * \param[in] tabafn Table potential file for angles
102 * \param[in] tabpfn Table potential file for proper dihedrals
103 * \param[in] tabbfn Table potential file for bonds
104 * \param[in] nbpu_opt Nonbonded Processing Unit (GPU/CPU etc.)
105 * \param[in] bNoSolvOpt Do not use solvent optimization
106 * \param[in] print_force Print forces for atoms with force >= print_force
108 void init_forcerec(FILE *fplog,
109 const output_env_t oenv,
110 t_forcerec *fr,
111 struct t_fcdata *fcd,
112 const t_inputrec *ir,
113 const gmx_mtop_t *mtop,
114 const t_commrec *cr,
115 matrix box,
116 const char *tabfn,
117 const char *tabafn,
118 const char *tabpfn,
119 const char *tabbfn,
120 const char *nbpu_opt,
121 gmx_bool bNoSolvOpt,
122 real print_force);
124 /*! \brief Divide exclusions over threads
126 * Set the exclusion load for the local exclusions and possibly threads
127 * \param[out] fr The force record
128 * \param[in] top The topology
130 void forcerec_set_excl_load(t_forcerec *fr,
131 const gmx_localtop_t *top);
133 /*! \brief Update parameters dependent on box
135 * Updates parameters in the forcerec that are time dependent
136 * \param[out] fr The force record
137 * \param[in] box The simulation box
139 void update_forcerec(t_forcerec *fr, matrix box);
141 #endif