Fixed include statements such that double precision version of genborn.c
[gromacs/rigid-bodies.git] / src / mdlib / wall.c
blob01a2a3005117de1c974e56c3b7841607681e797b
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
31 * And Hey:
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <math.h>
40 #include <string.h>
41 #include "maths.h"
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "macros.h"
45 #include "smalloc.h"
46 #include "force.h"
47 #include "filenm.h"
48 #include "nrnb.h"
49 #include "vec.h"
51 void make_wall_tables(FILE *fplog,const output_env_t oenv,
52 const t_inputrec *ir,const char *tabfn,
53 const gmx_groups_t *groups,
54 t_forcerec *fr)
56 int w,negp_pp,egp,i,j;
57 int *nm_ind;
58 char buf[STRLEN];
59 t_forcetable *tab;
61 negp_pp = ir->opts.ngener - ir->nwall;
62 nm_ind = groups->grps[egcENER].nm_ind;
64 if (fplog) {
65 fprintf(fplog,"Reading user tables for %d energy groups with %d walls\n",
66 negp_pp,ir->nwall);
69 snew(fr->wall_tab,ir->nwall);
70 for(w=0; w<ir->nwall; w++) {
71 snew(fr->wall_tab[w],negp_pp);
72 for(egp=0; egp<negp_pp; egp++) {
73 /* If the energy group pair is excluded, we don't need a table */
74 if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL)) {
75 tab = &fr->wall_tab[w][egp];
76 sprintf(buf,"%s",tabfn);
77 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,"_%s_%s.%s",
78 *groups->grpname[nm_ind[egp]],
79 *groups->grpname[nm_ind[negp_pp+w]],
80 ftp2ext(efXVG));
81 *tab = make_tables(fplog,oenv,fr,FALSE,buf,0,GMX_MAKETABLES_FORCEUSER);
82 /* Since wall have no charge, we can compress the table */
83 for(i=0; i<=tab->n; i++)
84 for(j=0; j<8; j++)
85 tab->tab[8*i+j] = tab->tab[12*i+4+j];
91 static void wall_error(int a,rvec *x,real r)
93 gmx_fatal(FARGS,
94 "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
95 "You might want to use the mdp option wall_r_linpot",
96 x[a][XX],x[a][YY],x[a][ZZ],r);
99 real do_walls(t_inputrec *ir,t_forcerec *fr,matrix box,t_mdatoms *md,
100 rvec x[],rvec f[],real lambda,real Vlj[],t_nrnb *nrnb)
102 int nwall,w,lam,i;
103 int ntw[2],at,ntype,ngid,ggid,*egp_flags,*type;
104 real *nbfp,lamfac,fac_d[2],fac_r[2],Cd,Cr,Vtot,Fwall[2];
105 real wall_z[2],r,mr,r1,r2,r4,Vd,Vr,V,Fd,Fr,F,dvdlambda;
106 dvec xf_z;
107 int n0,nnn;
108 real tabscale,*VFtab,rt,eps,eps2,Yt,Ft,Geps,Heps,Heps2,Fp,VV,FF;
109 unsigned short *gid=md->cENER;
110 t_forcetable *tab;
112 nwall = ir->nwall;
113 ngid = ir->opts.ngener;
114 ntype = fr->ntype;
115 nbfp = fr->nbfp;
116 egp_flags = fr->egp_flags;
118 for(w=0; w<nwall; w++) {
119 ntw[w] = 2*ntype*ir->wall_atomtype[w];
120 if (ir->wall_type == ewt93) {
121 fac_d[w] = ir->wall_density[w]*M_PI/6;
122 fac_r[w] = ir->wall_density[w]*M_PI/45;
123 } else if (ir->wall_type == ewt104) {
124 fac_d[w] = ir->wall_density[w]*M_PI/2;
125 fac_r[w] = ir->wall_density[w]*M_PI/5;
127 Fwall[w] = 0;
129 wall_z[0] = 0;
130 wall_z[1] = box[ZZ][ZZ];
132 Vtot = 0;
133 dvdlambda = 0;
134 clear_dvec(xf_z);
135 for(lam=0; lam<(md->nPerturbed ? 2 : 1); lam++) {
136 if (md->nPerturbed) {
137 if (lam == 0) {
138 lamfac = 1 - lambda;
139 type = md->typeA;
140 } else {
141 lamfac = 0;
142 type = md->typeB;
144 } else {
145 lamfac = 1;
146 type = md->typeA;
148 for(i=md->start; i<md->start+md->homenr; i++) {
149 for(w=0; w<nwall; w++) {
150 /* The wall energy groups are always at the end of the list */
151 ggid = gid[i]*ngid + ngid - nwall + w;
152 at = type[i];
153 Cd = nbfp[ntw[w]+2*at];
154 Cr = nbfp[ntw[w]+2*at+1];
155 if (!((Cd==0 && Cr==0) || egp_flags[ggid] & EGP_EXCL)) {
156 if (w == 0) {
157 r = x[i][ZZ];
158 } else {
159 r = wall_z[1] - x[i][ZZ];
161 if (r < ir->wall_r_linpot) {
162 mr = ir->wall_r_linpot - r;
163 r = ir->wall_r_linpot;
164 } else {
165 mr = 0;
167 if (ir->wall_type == ewtTABLE) {
168 if (r < 0)
169 wall_error(i,x,r);
170 tab = &(fr->wall_tab[w][gid[i]]);
171 tabscale = tab->scale;
172 VFtab = tab->tab;
174 rt = r*tabscale;
175 n0 = rt;
176 if (n0 >= tab->n) {
177 /* Beyond the table range, set V and F to zero */
178 V = 0;
179 F = 0;
180 } else {
181 eps = rt - n0;
182 eps2 = eps*eps;
183 /* Dispersion */
184 nnn = 8*n0;
185 Yt = VFtab[nnn];
186 Ft = VFtab[nnn+1];
187 Geps = VFtab[nnn+2]*eps;
188 Heps2 = VFtab[nnn+3]*eps2;
189 Fp = Ft + Geps + Heps2;
190 VV = Yt + Fp*eps;
191 FF = Fp + Geps + 2.0*Heps2;
192 Vd = Cd*VV;
193 Fd = Cd*FF;
194 /* Repulsion */
195 nnn = nnn + 4;
196 Yt = VFtab[nnn];
197 Ft = VFtab[nnn+1];
198 Geps = VFtab[nnn+2]*eps;
199 Heps2 = VFtab[nnn+3]*eps2;
200 Fp = Ft + Geps + Heps2;
201 VV = Yt + Fp*eps;
202 FF = Fp + Geps + 2.0*Heps2;
203 Vr = Cr*VV;
204 Fr = Cr*FF;
205 V = Vd + Vr;
206 F = -lamfac*(Fd + Fr)*tabscale;
208 } else if (ir->wall_type == ewt93) {
209 if (r <= 0)
210 wall_error(i,x,r);
211 r1 = 1/r;
212 r2 = r1*r1;
213 r4 = r2*r2;
214 Vd = fac_d[w]*Cd*r2*r1;
215 Vr = fac_r[w]*Cr*r4*r4*r1;
216 V = Vr - Vd;
217 F = lamfac*(9*Vr - 3*Vd)*r1;
218 } else {
219 if (r <= 0)
220 wall_error(i,x,r);
221 r1 = 1/r;
222 r2 = r1*r1;
223 r4 = r2*r2;
224 Vd = fac_d[w]*Cd*r4;
225 Vr = fac_r[w]*Cr*r4*r4*r2;
226 V = Vr - Vd;
227 F = lamfac*(10*Vr - 4*Vd)*r1;
229 if (mr > 0) {
230 V += mr*F;
232 if (w == 1) {
233 F = -F;
235 Vlj[ggid] += lamfac*V;
236 Vtot += V;
237 f[i][ZZ] += F;
238 /* Because of the single sum virial calculation we need to add
239 * the full virial contribution of the walls.
240 * Since the force only has a z-component, there is only
241 * a contribution to the z component of the virial tensor.
242 * We could also determine the virial contribution directly,
243 * which would be cheaper here, but that would require extra
244 * communication for f_novirsum for with virtual sites in parallel.
246 xf_z[XX] -= x[i][XX]*F;
247 xf_z[YY] -= x[i][YY]*F;
248 xf_z[ZZ] -= wall_z[w]*F;
252 if (md->nPerturbed)
253 dvdlambda += (lam==0 ? -1 : 1)*Vtot;
255 inc_nrnb(nrnb,eNR_WALLS,md->homenr);
258 for(i=0; i<DIM; i++) {
259 fr->vir_wall_z[i] = -0.5*xf_z[i];
262 return dvdlambda;