Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / mdlib / adress.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
5  * Copyright (c) 2011,2014, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief Implementation of the AdResS method.
38  */
39 #ifndef _adress_h_
40 #define _adress_h_
41
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44
45 #ifdef __cplusplus
46 extern "C"
47 {
48 #endif
49
50 struct t_pbc;
51
52 /** \brief calculates the AdResS weight of a particle
53  *
54  * \param[in] x position of the particle
55  * \param[in] adresstype type of address weight function
56  *                       eAdressOff - explicit simulation
57  *                       eAdressConst - constant weight all over the box
58  *                       eAdressXSplit - split in x direction with ref as center
59  *                       eAdressSphere - spherical splitting with ref as center
60  *                       else - weight = 1 - explicit simulation
61  * \param[in] adressr radius/size of the explicit zone
62  * \param[in] adressw size of the hybrid zone
63  * \param[in] ref center of the explicit zone
64  *                for adresstype 1 - unused
65  *                for adresstype 2 - only ref[0] is used
66  * \param[in] pbc pbc struct for calculating shortest distance
67  * \param[in] fr the forcerec containing all the parameters
68  *
69  * \return weight of the particle
70  *
71  */
72 real
73 adress_weight(rvec                 x,
74               int                  adresstype,
75               real                 adressr,
76               real                 adressw,
77               rvec     *           ref,
78               struct t_pbc     *   pbc,
79               t_forcerec *         fr);
80
81 /** \brief update the weight of all coarse-grained particles in several charge groups for com vsites
82  *
83  * \param[in,out] fplog log file in case of debug
84  * \param[in] cg0 first charge group to update
85  * \param[in] cg1 last+1 charge group to update
86  * \param[in] cgs block containing the cg index
87  * \param[in] x array with all the particle positions
88  * \param[in] fr the forcerec containing all the parameters
89  * \param[in,out] mdatoms the struct containing all the atoms properties
90  * \param[in] pbc for shortest distance in adress_weight
91  */
92 void
93 update_adress_weights_com(FILE *               fplog,
94                           int                  cg0,
95                           int                  cg1,
96                           t_block *            cgs,
97                           rvec                 x[],
98                           t_forcerec *         fr,
99                           t_mdatoms *          mdatoms,
100                           struct t_pbc *       pbc);
101
102 /** \brief update the weight of all coarse-grained particles for cog vsites
103  *
104  * \param[in] ip contains interaction parameters, in this case the number of constructing atoms n for vsitesn
105  * \param[in] ilist list of interaction types, in this case the virtual site types are what's important
106  * \param[in] x array with all the particle positions
107  * \param[in] fr the forcerec containing all the parameters
108  * \param[in,out] mdatoms the struct containing all the atoms properties
109  * \param[in] pbc for shortest distance in adress_weight
110  */
111 void
112 update_adress_weights_cog(t_iparams            ip[],
113                           t_ilist              ilist[],
114                           rvec                 x[],
115                           t_forcerec *         fr,
116                           t_mdatoms *          mdatoms,
117                           struct t_pbc *       pbc);
118
119 /** \brief update the weight of all coarse-grained particles in several charge groups for atom vsites
120  *
121  * \param[in] cg0 first charge group to update
122  * \param[in] cg1 last+1 charge group to update
123  * \param[in] cgs block containing the cg index
124  * \param[in] x array with all the particle positions
125  * \param[in] fr the forcerec containing all the parameters
126  * \param[in,out] mdatoms the struct containing all the atoms properties
127  * \param[in] pbc for shortest distance in adress_weight
128  */
129 void
130 update_adress_weights_atom(int                  cg0,
131                            int                  cg1,
132                            t_block *            cgs,
133                            rvec                 x[],
134                            t_forcerec *         fr,
135                            t_mdatoms *          mdatoms,
136                            struct t_pbc *       pbc);
137
138 /** \brief update the weight on per atom basis of all coarse-grained particles in several charge groups for atom vsites
139  *
140  * \param[in] cg0 first charge group to update
141  * \param[in] cg1 last+1 charge group to update
142  * \param[in] cgs block containing the cg index
143  * \param[in] x array with all the particle positions
144  * \param[in] fr the forcerec containing all the parameters
145  * \param[in,out] mdatoms the struct containing all the atoms properties
146  * \param[in] pbc for shortest distance in adress_weight
147  */
148 void
149 update_adress_weights_atom_per_atom(int                  cg0,
150                                     int                  cg1,
151                                     t_block *            cgs,
152                                     rvec                 x[],
153                                     t_forcerec *         fr,
154                                     t_mdatoms *          mdatoms,
155                                     struct t_pbc *       pbc);
156
157 /** \brief add AdResS IC thermodynamic force to f_novirsum
158  *
159  * \param[in] cg0 first charge group to update
160  * \param[in] cg1 last+1 charge group to update
161  * \param[in] cgs block containing the cg index
162  * \param[in] x array with all the particle positions
163  * \param[in,out] f the force array pointing at f_novirsum from sim_util.c
164  * \param[in] fr the forcerec containing all the parameters
165  * \param[in] mdatoms the struct containing all the atoms properties
166  * \param[in] pbc for shortest distance to fr->adress_refs
167  */
168 void
169 adress_thermo_force(int                  cg0,
170                     int                  cg1,
171                     t_block *            cgs,
172                     rvec                 x[],
173                     rvec                 f[],
174                     t_forcerec *         fr,
175                     t_mdatoms *          mdatoms,
176                     struct t_pbc *       pbc);
177
178
179 /** \brief checks weather a cpu calculates only coarse-grained or explicit interactions
180  *
181  * \param[in] n_ex number of explicit particles
182  * \param[in] n_hyb number of hybrid particles
183  * \param[in] n_cg number of coarse-grained particles
184  * \param[in,out] mdatoms the struct containing all the atoms properties
185  */
186 void adress_set_kernel_flags(int n_ex, int n_hyb, int n_cg, t_mdatoms * mdatoms);
187
188 /** \brief looks up  if a energy group is explicit
189  * \param[in] fr the forcerec containing all the parameters
190  * \param[in] egp_nr energy group number
191  * \return boolean if explicit or not
192  */
193 gmx_bool egp_explicit(t_forcerec *   fr, int egp_nr);
194
195 /** \brief looks up  if a energy group is coarse-grained
196  * \param[in] fr the forcerec containing all the parameters
197  * \param[in] egp_nr energy group number
198  * \return boolean if coarse-grained or not
199  */
200 gmx_bool egp_coarsegrained(t_forcerec *   fr, int egp_nr);
201
202 #ifdef __cplusplus
203 }
204 #endif
205
206 #endif