Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / mdlib / adress.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 4.0.5
11  * Written by Christoph Junghans, Brad Lambeth, and possibly others.
12  * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
13  * All rights reserved.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35
36 #ifndef _adress_h_
37 #define _adress_h_
38
39 /** \file adress.h
40  *
41  * \brief Implementation of the AdResS method
42  *
43  */
44
45 #include "types/simple.h"
46 #include "typedefs.h"
47
48 /** \brief calculates the AdResS weight of a particle
49  *
50  * \param[in] x position of the particle
51  * \param[in] adresstype type of address weight function
52  *                       eAdressOff - explicit simulation
53  *                       eAdressConst - constant weight all over the box
54  *                       eAdressXSplit - split in x direction with ref as center
55  *                       eAdressSphere - spherical splitting with ref as center
56  *                       else - weight = 1 - explicit simulation
57  * \param[in] adressr radius/size of the explicit zone
58  * \param[in] adressw size of the hybrid zone
59  * \param[in] ref center of the explicit zone
60  *                for adresstype 1 - unused
61  *                for adresstype 2 - only ref[0] is used
62  * \param[in] pbc pbc struct for calculating shortest distance
63  * \param[in] fr the forcerec containing all the parameters
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
115 /** \brief update the weight of all coarse-grained particles in several charge groups for atom vsites
116  *
117  * \param[in] cg0 first charge group to update
118  * \param[in] cg1 last+1 charge group to update
119  * \param[in] cgs block containing the cg index
120  * \param[in] x array with all the particle positions
121  * \param[in] fr the forcerec containing all the parameters
122  * \param[in,out] mdatoms the struct containing all the atoms properties
123  * \param[in] pbc for shortest distance in adress_weight
124  */
125 void
126 update_adress_weights_atom(int                  cg0,
127                            int                  cg1,
128                            t_block *            cgs,
129                            rvec                 x[],
130                            t_forcerec *         fr,
131                            t_mdatoms *          mdatoms,
132                            t_pbc *              pbc);
133
134 /** \brief update the weight on per atom basis of all coarse-grained particles in several charge groups for atom vsites
135  *
136  * \param[in] cg0 first charge group to update
137  * \param[in] cg1 last+1 charge group to update
138  * \param[in] cgs block containing the cg index
139  * \param[in] x array with all the particle positions
140  * \param[in] fr the forcerec containing all the parameters
141  * \param[in,out] mdatoms the struct containing all the atoms properties
142  * \param[in] pbc for shortest distance in adress_weight
143  */
144 void
145 update_adress_weights_atom_per_atom(int                  cg0,
146                                     int                  cg1,
147                                     t_block *            cgs,
148                                     rvec                 x[],
149                                     t_forcerec *         fr,
150                                     t_mdatoms *          mdatoms,
151                                     t_pbc *              pbc);
152
153 /** \brief add AdResS IC thermodynamic force to f_novirsum
154  *
155  * \param[in] cg0 first charge group to update
156  * \param[in] cg1 last+1 charge group to update
157  * \param[in] cgs block containing the cg index
158  * \param[in] x array with all the particle positions
159  * \param[in,out] f the force array pointing at f_novirsum from sim_util.c
160  * \param[in] fr the forcerec containing all the parameters
161  * \param[in] mdatoms the struct containing all the atoms properties
162  * \param[in] pbc for shortest distance to fr->adress_refs
163  */
164 void
165 adress_thermo_force(int                  cg0,
166                     int                  cg1,
167                     t_block *            cgs,
168                     rvec                 x[],
169                     rvec                 f[],
170                     t_forcerec *         fr,
171                     t_mdatoms *          mdatoms,
172                     t_pbc *              pbc);
173
174
175 /** \brief checks weather a cpu calculates only coarse-grained or explicit interactions
176  *
177  * \param[in] n_ex number of explicit particles
178  * \param[in] n_hyb number of hybrid particles
179  * \param[in] n_cg number of coarse-grained particles
180  * \param[in,out] mdatoms the struct containing all the atoms properties
181  */
182 void adress_set_kernel_flags(int n_ex, int n_hyb, int n_cg, t_mdatoms * mdatoms);
183
184 /** \brief looks up  if a energy group is explicit
185  * \param[in] fr the forcerec containing all the parameters
186  * \param[in] egp_nr energy group number
187  * \return boolean if explicit or not
188  */
189 gmx_bool egp_explicit(t_forcerec *   fr, int egp_nr);
190
191 /** \brief looks up  if a energy group is coarse-grained
192  * \param[in] fr the forcerec containing all the parameters
193  * \param[in] egp_nr energy group number
194  * \return boolean if coarse-grained or not
195  */
196 gmx_bool egp_coarsegrained(t_forcerec *   fr, int egp_nr);
197 #endif