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