Merge release-4-6 into master
[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, 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
37 #ifndef _adress_h_
38 #define _adress_h_
39
40 /** \file adress.h
41  *
42  * \brief Implementation of the AdResS method
43  *
44  */
45
46 #include "types/simple.h"
47 #include "typedefs.h"
48
49 /** \brief calculates the AdResS weight of a particle
50  *
51  * \param[in] x position of the particle
52  * \param[in] adresstype type of address weight function
53  *                       eAdressOff - explicit simulation
54  *                       eAdressConst - constant weight all over the box
55  *                       eAdressXSplit - split in x direction with ref as center
56  *                       eAdressSphere - spherical splitting with ref as center
57  *                       else - weight = 1 - explicit simulation
58  * \param[in] adressr radius/size of the explicit zone
59  * \param[in] adressw size of the hybrid zone
60  * \param[in] ref center of the explicit zone
61  *                for adresstype 1 - unused
62  *                for adresstype 2 - only ref[0] is used
63  * \param[in] pbc pbc struct for calculating shortest distance
64  * \param[in] fr the forcerec containing all the parameters
65  *
66  * \return weight of the particle
67  *
68  */
69 real
70 adress_weight(rvec                 x,
71               int                  adresstype,
72               real                 adressr,
73               real                 adressw,
74               rvec     *           ref,
75               t_pbc     *          pbc,
76               t_forcerec *         fr);
77
78 /** \brief update the weight of all coarse-grained particles in several charge groups for com vsites
79  *
80  * \param[in,out] fplog log file in case of debug
81  * \param[in] cg0 first charge group to update
82  * \param[in] cg1 last+1 charge group to update
83  * \param[in] cgs block containing the cg index
84  * \param[in] x array with all the particle positions
85  * \param[in] fr the forcerec containing all the parameters
86  * \param[in,out] mdatoms the struct containing all the atoms properties
87  * \param[in] pbc for shortest distance in adress_weight
88  */
89 void
90 update_adress_weights_com(FILE *               fplog,
91                           int                  cg0,
92                           int                  cg1,
93                           t_block *            cgs,
94                           rvec                 x[],
95                           t_forcerec *         fr,
96                           t_mdatoms *          mdatoms,
97                           t_pbc *              pbc);
98
99 /** \brief update the weight of all coarse-grained particles for cog vsites
100  *
101  * \param[in] ip contains interaction parameters, in this case the number of constructing atoms n for vsitesn
102  * \param[in] ilist list of interaction types, in this case the virtual site types are what's important
103  * \param[in] x array with all the particle positions
104  * \param[in] fr the forcerec containing all the parameters
105  * \param[in,out] mdatoms the struct containing all the atoms properties
106  * \param[in] pbc for shortest distance in adress_weight
107  */
108 void
109 update_adress_weights_cog(t_iparams            ip[],
110                           t_ilist              ilist[],
111                           rvec                 x[],
112                           t_forcerec *         fr,
113                           t_mdatoms *          mdatoms,
114                           t_pbc *              pbc);
115
116 /** \brief update the weight of all coarse-grained particles in several charge groups for atom vsites
117  *
118  * \param[in] cg0 first charge group to update
119  * \param[in] cg1 last+1 charge group to update
120  * \param[in] cgs block containing the cg index
121  * \param[in] x array with all the particle positions
122  * \param[in] fr the forcerec containing all the parameters
123  * \param[in,out] mdatoms the struct containing all the atoms properties
124  * \param[in] pbc for shortest distance in adress_weight
125  */
126 void
127 update_adress_weights_atom(int                  cg0,
128                            int                  cg1,
129                            t_block *            cgs,
130                            rvec                 x[],
131                            t_forcerec *         fr,
132                            t_mdatoms *          mdatoms,
133                            t_pbc *              pbc);
134
135 /** \brief update the weight on per atom basis of all coarse-grained particles in several charge groups for atom vsites
136  *
137  * \param[in] cg0 first charge group to update
138  * \param[in] cg1 last+1 charge group to update
139  * \param[in] cgs block containing the cg index
140  * \param[in] x array with all the particle positions
141  * \param[in] fr the forcerec containing all the parameters
142  * \param[in,out] mdatoms the struct containing all the atoms properties
143  * \param[in] pbc for shortest distance in adress_weight
144  */
145 void
146 update_adress_weights_atom_per_atom(int                  cg0,
147                                     int                  cg1,
148                                     t_block *            cgs,
149                                     rvec                 x[],
150                                     t_forcerec *         fr,
151                                     t_mdatoms *          mdatoms,
152                                     t_pbc *              pbc);
153
154 /** \brief add AdResS IC thermodynamic force to f_novirsum
155  *
156  * \param[in] cg0 first charge group to update
157  * \param[in] cg1 last+1 charge group to update
158  * \param[in] cgs block containing the cg index
159  * \param[in] x array with all the particle positions
160  * \param[in,out] f the force array pointing at f_novirsum from sim_util.c
161  * \param[in] fr the forcerec containing all the parameters
162  * \param[in] mdatoms the struct containing all the atoms properties
163  * \param[in] pbc for shortest distance to fr->adress_refs
164  */
165 void
166 adress_thermo_force(int                  cg0,
167                     int                  cg1,
168                     t_block *            cgs,
169                     rvec                 x[],
170                     rvec                 f[],
171                     t_forcerec *         fr,
172                     t_mdatoms *          mdatoms,
173                     t_pbc *              pbc);
174
175
176 /** \brief checks weather a cpu calculates only coarse-grained or explicit interactions
177  *
178  * \param[in] n_ex number of explicit particles
179  * \param[in] n_hyb number of hybrid particles
180  * \param[in] n_cg number of coarse-grained particles
181  * \param[in,out] mdatoms the struct containing all the atoms properties
182  */
183 void adress_set_kernel_flags(int n_ex, int n_hyb, int n_cg, t_mdatoms * mdatoms);
184
185 /** \brief looks up  if a energy group is explicit
186  * \param[in] fr the forcerec containing all the parameters
187  * \param[in] egp_nr energy group number
188  * \return boolean if explicit or not
189  */
190 gmx_bool egp_explicit(t_forcerec *   fr, int egp_nr);
191
192 /** \brief looks up  if a energy group is coarse-grained
193  * \param[in] fr the forcerec containing all the parameters
194  * \param[in] egp_nr energy group number
195  * \return boolean if coarse-grained or not
196  */
197 gmx_bool egp_coarsegrained(t_forcerec *   fr, int egp_nr);
198 #endif