Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / 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) 2012,2013, by the GROMACS development team, led by
6  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
7  * others, as listed in the AUTHORS file in the top-level source
8  * 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 for calculating shortest distance to ref
64  *
65  * \return weight of the particle
66  *
67  */
68 real 
69 adress_weight(rvec             x,
70               int              adresstype,
71               real             adressr,
72               real             adressw,
73               rvec *           ref,
74               t_pbc *          pbc,
75               t_forcerec *         fr);
76
77 /** \brief update the weight of all coarse-grained particles in several charge groups for com vsites
78  *
79  * \param[in,out] fplog log file in case of debug
80  * \param[in] cg0 first charge group to update
81  * \param[in] cg1 last+1 charge group to update
82  * \param[in] cgs block containing the cg index 
83  * \param[in] x array with all the particle positions  
84  * \param[in] fr the forcerec containing all the parameters
85  * \param[in,out] mdatoms the struct containing all the atoms properties
86  * \param[in] pbc for shortest distance in adress_weight
87  */
88 void
89 update_adress_weights_com(FILE *               fplog,
90                           int                  cg0,
91                           int                  cg1,
92                           t_block *            cgs,
93                           rvec                 x[],
94                           t_forcerec *         fr,
95                           t_mdatoms *          mdatoms,
96                           t_pbc *              pbc);
97
98 /** \brief update the weight of all coarse-grained particles for cog vsites
99  *
100  * \param[in] ip contains interaction parameters, in this case the number of constructing atoms n for vsitesn
101  * \param[in] ilist list of interaction types, in this case the virtual site types are what's important
102  * \param[in] x array with all the particle positions  
103  * \param[in] fr the forcerec containing all the parameters
104  * \param[in,out] mdatoms the struct containing all the atoms properties
105  * \param[in] pbc for shortest distance in adress_weight
106  */
107 void
108 update_adress_weights_cog(t_iparams            ip[],
109                           t_ilist              ilist[],
110                           rvec                 x[],
111                           t_forcerec *         fr,
112                           t_mdatoms *          mdatoms,
113                           t_pbc *              pbc);
114 /** \brief update the weight of all coarse-grained particles in several charge groups for atom vsites
115  *
116  * \param[in] cg0 first charge group to update
117  * \param[in] cg1 last+1 charge group to update
118  * \param[in] cgs block containing the cg index 
119  * \param[in] x array with all the particle positions  
120  * \param[in] fr the forcerec containing all the parameters
121  * \param[in,out] mdatoms the struct containing all the atoms properties
122  * \param[in] pbc for shortest distance in adress_weight
123  */
124 void
125 update_adress_weights_atom(int                  cg0,
126                            int                  cg1,
127                            t_block *            cgs,
128                            rvec                 x[],
129                            t_forcerec *         fr,
130                            t_mdatoms *          mdatoms,
131                            t_pbc *              pbc);
132
133 void
134 update_adress_weights_atom_per_atom(int                  cg0,
135                            int                  cg1,
136                            t_block *            cgs,
137                            rvec                 x[],
138                            t_forcerec *         fr,
139                            t_mdatoms *          mdatoms,
140                            t_pbc *              pbc);
141
142 /** \brief add AdResS IC thermodynamic force to f_novirsum
143  *
144  * \param[in] cg0 first charge group to update
145  * \param[in] cg1 last+1 charge group to update
146  * \param[in] cgs block containing the cg index 
147  * \param[in] x array with all the particle positions  
148  * \param[in,out] f the force array pointing at f_novirsum from sim_util.c
149  * \param[in] fr the forcerec containing all the parameters
150  * \param[in] mdatoms the struct containing all the atoms properties
151  * \param[in] pbc for shortest distance to fr->adress_refs
152  */
153 void
154 adress_thermo_force(int                  cg0,
155                     int                  cg1,
156                     t_block *            cgs,
157                     rvec                 x[],
158                     rvec                 f[],
159                     t_forcerec *         fr,
160                     t_mdatoms *          mdatoms,
161                     t_pbc *              pbc);
162
163 void adress_set_kernel_flags(int n_ex, int n_hyb, int n_cg, t_mdatoms * mdatoms);
164
165 /* functions to look up if a energy group is explicit or coarse-grained*/
166 gmx_bool egp_explicit(t_forcerec *   fr, int egp_nr);
167 gmx_bool egp_coarsegrained(t_forcerec *   fr, int egp_nr);
168 #endif