Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / mdlib / adress.h
blobc14fd10dac03d1ffb7f40992d26a339e373021d7
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifndef _adress_h_
37 #define _adress_h_
39 /** \file adress.h
41 * \brief Implementation of the AdResS method
45 #include "types/simple.h"
46 #include "typedefs.h"
48 /** \brief calculates the AdResS weight of a particle
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 for calculating shortest distance to ref
64 * \return weight of the particle
67 real
68 adress_weight(rvec x,
69 int adresstype,
70 real adressr,
71 real adressw,
72 rvec * ref,
73 t_pbc * pbc,
74 t_forcerec * fr);
76 /** \brief update the weight of all coarse-grained particles in several charge groups for com vsites
78 * \param[in,out] fplog log file in case of debug
79 * \param[in] cg0 first charge group to update
80 * \param[in] cg1 last+1 charge group to update
81 * \param[in] cgs block containing the cg index
82 * \param[in] x array with all the particle positions
83 * \param[in] fr the forcerec containing all the parameters
84 * \param[in,out] mdatoms the struct containing all the atoms properties
85 * \param[in] pbc for shortest distance in adress_weight
87 void
88 update_adress_weights_com(FILE * fplog,
89 int cg0,
90 int cg1,
91 t_block * cgs,
92 rvec x[],
93 t_forcerec * fr,
94 t_mdatoms * mdatoms,
95 t_pbc * pbc);
97 /** \brief update the weight of all coarse-grained particles for cog vsites
99 * \param[in] ip contains interaction parameters, in this case the number of constructing atoms n for vsitesn
100 * \param[in] ilist list of interaction types, in this case the virtual site types are what's important
101 * \param[in] x array with all the particle positions
102 * \param[in] fr the forcerec containing all the parameters
103 * \param[in,out] mdatoms the struct containing all the atoms properties
104 * \param[in] pbc for shortest distance in adress_weight
106 void
107 update_adress_weights_cog(t_iparams ip[],
108 t_ilist ilist[],
109 rvec x[],
110 t_forcerec * fr,
111 t_mdatoms * mdatoms,
112 t_pbc * pbc);
113 /** \brief update the weight of all coarse-grained particles in several charge groups for atom vsites
115 * \param[in] cg0 first charge group to update
116 * \param[in] cg1 last+1 charge group to update
117 * \param[in] cgs block containing the cg index
118 * \param[in] x array with all the particle positions
119 * \param[in] fr the forcerec containing all the parameters
120 * \param[in,out] mdatoms the struct containing all the atoms properties
121 * \param[in] pbc for shortest distance in adress_weight
123 void
124 update_adress_weights_atom(int cg0,
125 int cg1,
126 t_block * cgs,
127 rvec x[],
128 t_forcerec * fr,
129 t_mdatoms * mdatoms,
130 t_pbc * pbc);
132 void
133 update_adress_weights_atom_per_atom(int cg0,
134 int cg1,
135 t_block * cgs,
136 rvec x[],
137 t_forcerec * fr,
138 t_mdatoms * mdatoms,
139 t_pbc * pbc);
141 /** \brief add AdResS IC thermodynamic force to f_novirsum
143 * \param[in] cg0 first charge group to update
144 * \param[in] cg1 last+1 charge group to update
145 * \param[in] cgs block containing the cg index
146 * \param[in] x array with all the particle positions
147 * \param[in,out] f the force array pointing at f_novirsum from sim_util.c
148 * \param[in] fr the forcerec containing all the parameters
149 * \param[in] mdatoms the struct containing all the atoms properties
150 * \param[in] pbc for shortest distance to fr->adress_refs
152 void
153 adress_thermo_force(int cg0,
154 int cg1,
155 t_block * cgs,
156 rvec x[],
157 rvec f[],
158 t_forcerec * fr,
159 t_mdatoms * mdatoms,
160 t_pbc * pbc);
162 void adress_set_kernel_flags(int n_ex, int n_hyb, int n_cg, t_mdatoms * mdatoms);
164 /* functions to look up if a energy group is explicit or coarse-grained*/
165 gmx_bool egp_explicit(t_forcerec * fr, int egp_nr);
166 gmx_bool egp_coarsegrained(t_forcerec * fr, int egp_nr);
167 #endif