3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
58 /****************************************************************************/
59 /* This program calculates the ordering of water molecules across a box, as */
60 /* function of the z-coordinate. This implies averaging over slices and over*/
61 /* time. Output is the average cos of the angle of the dipole with the */
62 /* normal, per slice. */
63 /* In addition, it calculates the average dipole moment itself in three */
65 /****************************************************************************/
67 void calc_h2order(const char *fn
, atom_id index
[], int ngx
, rvec
**slDipole
,
68 real
**slOrder
, real
*slWidth
, int *nslices
,
69 t_topology
*top
, int ePBC
,
70 int axis
, bool bMicel
, atom_id micel
[], int nmic
,
71 const output_env_t oenv
)
73 rvec
*x0
, /* coordinates with pbc */
74 dipole
, /* dipole moment due to one molecules */
76 com
; /* center of mass of micel, with bMicel */
77 rvec
*dip
; /* sum of dipoles, unnormalized */
78 matrix box
; /* box (3x3) */
80 real t
, /* time from trajectory */
81 *sum
, /* sum of all cosines of dipoles, per slice */
82 *frame
; /* order over one frame */
83 int natoms
, /* nr. atoms in trj */
85 slice
=0, /* current slice number */
86 *count
; /* nr. of atoms in one slice */
88 if ((natoms
= read_first_x(oenv
,&status
,fn
,&t
,&x0
,box
)) == 0)
89 gmx_fatal(FARGS
,"Could not read coordinates from statusfile\n");
92 *nslices
= (int)(box
[axis
][axis
] * 10); /* default value */
98 normal
[0] = 1; normal
[1] = 0; normal
[2] = 0;
101 normal
[0] = 0; normal
[1] = 1; normal
[2] = 0;
104 normal
[0] = 0; normal
[1] = 0; normal
[2] = 1;
107 gmx_fatal(FARGS
,"No valid value for -axis-. Exiting.\n");
108 /* make compiler happy */
109 normal
[0] = 1; normal
[1] = 0; normal
[2] = 0;
113 snew(count
, *nslices
);
116 snew(frame
, *nslices
);
118 *slWidth
= box
[axis
][axis
]/(*nslices
);
119 fprintf(stderr
,"Box divided in %d slices. Initial width of slice: %f\n",
124 /*********** Start processing trajectory ***********/
127 *slWidth
= box
[axis
][axis
]/(*nslices
);
130 rm_pbc(&(top
->idef
),ePBC
,top
->atoms
.nr
,box
,x0
,x0
);
133 calc_xcm(x0
, nmic
, micel
, top
->atoms
.atom
, com
, FALSE
);
135 for (i
= 0; i
< ngx
/3; i
++)
137 /* put all waters in box */
138 for (j
= 0; j
< DIM
; j
++)
140 if (x0
[index
[3*i
]][j
] < 0)
142 x0
[index
[3*i
]][j
] += box
[j
][j
];
143 x0
[index
[3*i
+1]][j
] += box
[j
][j
];
144 x0
[index
[3*i
+2]][j
] += box
[j
][j
];
146 if (x0
[index
[3*i
]][j
] > box
[j
][j
])
148 x0
[index
[3*i
]][j
] -= box
[j
][j
];
149 x0
[index
[3*i
+1]][j
] -= box
[j
][j
];
150 x0
[index
[3*i
+2]][j
] -= box
[j
][j
];
154 for (j
= 0; j
< DIM
; j
++)
156 x0
[index
[3*i
]][j
] * top
->atoms
.atom
[index
[3*i
]].q
+
157 x0
[index
[3*i
+1]][j
] * top
->atoms
.atom
[index
[3*i
+1]].q
+
158 x0
[index
[3*i
+2]][j
] * top
->atoms
.atom
[index
[3*i
+2]].q
;
160 /* now we have a dipole vector. Might as well safe it. Then the
161 rest depends on whether we're dealing with
162 a flat or a spherical interface.
166 { /* this is for spherical interfaces */
167 rvec_sub(com
, x0
[index
[3*i
]], normal
); /* vector from Oxygen to COM */
168 slice
= norm(normal
)/(*slWidth
); /* spherical slice */
170 sum
[slice
] += iprod(dipole
, normal
) / (norm(dipole
) * norm(normal
));
171 frame
[slice
] += iprod(dipole
, normal
) / (norm(dipole
) * norm(normal
));
175 /* this is for flat interfaces */
177 /* determine which slice atom is in */
178 slice
= (x0
[index
[3*i
]][axis
] / (*slWidth
));
179 if (slice
< 0 || slice
>= *nslices
)
181 fprintf(stderr
,"Coordinate: %f ",x0
[index
[3*i
]][axis
]);
182 fprintf(stderr
,"HELP PANIC! slice = %d, OUT OF RANGE!\n",slice
);
186 rvec_add(dipole
, dip
[slice
], dip
[slice
]);
187 /* Add dipole to total. mag[slice] is total dipole in axis direction */
188 sum
[slice
] += iprod(dipole
,normal
)/norm(dipole
);
189 frame
[slice
] += iprod(dipole
,normal
)/norm(dipole
);
190 /* increase count for that slice */
196 } while (read_next_x(oenv
,status
,&t
,natoms
,x0
,box
));
197 /*********** done with status file **********/
199 fprintf(stderr
,"\nRead trajectory. Printing parameters to file\n");
201 for (i
= 0; i
< *nslices
; i
++) /* average over frames */
203 fprintf(stderr
,"%d waters in slice %d\n",count
[i
],i
);
204 if (count
[i
] > 0) /* divide by number of molecules in each slice */
206 sum
[i
] = sum
[i
] / count
[i
];
207 dip
[i
][XX
] = dip
[i
][XX
] / count
[i
];
208 dip
[i
][YY
] = dip
[i
][YY
] / count
[i
];
209 dip
[i
][ZZ
] = dip
[i
][ZZ
] / count
[i
];
212 fprintf(stderr
,"No water in slice %d\n",i
);
215 *slOrder
= sum
; /* copy a pointer, I hope */
217 sfree(x0
); /* free memory used by coordinate arrays */
220 void h2order_plot(rvec dipole
[], real order
[], const char *afile
,
221 int nslices
, real slWidth
, const output_env_t oenv
)
223 FILE *ord
; /* xvgr files with order parameters */
224 int slice
; /* loop index */
225 char buf
[256]; /* for xvgr title */
226 real factor
; /* conversion to Debye from electron*nm */
228 /* factor = 1e-9*1.60217733e-19/3.336e-30 */
229 factor
= 1.60217733/3.336e-2;
230 fprintf(stderr
,"%d slices\n",nslices
);
231 sprintf(buf
,"Water orientation with respect to normal");
232 ord
= xvgropen(afile
,buf
,
233 "box (nm)","mu_x, mu_y, mu_z (D), cosine with normal",oenv
);
235 for (slice
= 0; slice
< nslices
; slice
++)
236 fprintf(ord
,"%8.3f %8.3f %8.3f %8.3f %e\n", slWidth
*slice
,
237 factor
*dipole
[slice
][XX
], factor
*dipole
[slice
][YY
],
238 factor
*dipole
[slice
][ZZ
], order
[slice
]);
243 int gmx_h2order(int argc
,char *argv
[])
245 const char *desc
[] = {
246 "Compute the orientation of water molecules with respect to the normal",
247 "of the box. The program determines the average cosine of the angle",
248 "between de dipole moment of water and an axis of the box. The box is",
249 "divided in slices and the average orientation per slice is printed.",
250 "Each water molecule is assigned to a slice, per time frame, based on the",
251 "position of the oxygen. When -nm is used the angle between the water",
252 "dipole and the axis from the center of mass to the oxygen is calculated",
253 "instead of the angle between the dipole and a box axis."
255 static int axis
= 2; /* normal to memb. default z */
256 static const char *axtitle
="Z";
257 static int nslices
= 0; /* nr of slices defined */
259 { "-d", FALSE
, etSTR
, {&axtitle
},
260 "Take the normal on the membrane in direction X, Y or Z." },
261 { "-sl", FALSE
, etINT
, {&nslices
},
262 "Calculate order parameter as function of boxlength, dividing the box"
265 const char *bugs
[] = {
266 "The program assigns whole water molecules to a slice, based on the first"
267 "atom of three in the index file group. It assumes an order O,H,H."
268 "Name is not important, but the order is. If this demand is not met,"
269 "assigning molecules to slices is different."
273 real
*slOrder
, /* av. cosine, per slice */
274 slWidth
= 0.0; /* width of a slice */
276 char *grpname
, /* groupnames */
278 int ngx
, /* nr. of atomsin sol group */
279 nmic
; /* nr. of atoms in micelle */
280 t_topology
*top
; /* topology */
282 atom_id
*index
, /* indices for solvent group */
284 bool bMicel
= FALSE
; /* think we're a micel */
285 t_filenm fnm
[] = { /* files for g_order */
286 { efTRX
, "-f", NULL
, ffREAD
}, /* trajectory file */
287 { efNDX
, NULL
, NULL
, ffREAD
}, /* index file */
288 { efNDX
, "-nm", NULL
, ffOPTRD
}, /* index with micelle atoms */
289 { efTPX
, NULL
, NULL
, ffREAD
}, /* topology file */
290 { efXVG
,"-o", "order", ffWRITE
}, /* xvgr output file */
293 #define NFILE asize(fnm)
295 CopyRight(stderr
,argv
[0]);
296 parse_common_args(&argc
, argv
,
297 PCA_CAN_VIEW
| PCA_CAN_TIME
| PCA_BE_NICE
, NFILE
,
298 fnm
, asize(pa
),pa
,asize(desc
),desc
,asize(bugs
),bugs
,&oenv
);
299 bMicel
= opt2bSet("-nm",NFILE
,fnm
);
301 top
= read_top(ftp2fn(efTPX
,NFILE
,fnm
),&ePBC
); /* read topology file */
303 rd_index(ftp2fn(efNDX
,NFILE
,fnm
),1,&ngx
,&index
,&grpname
);
306 rd_index(opt2fn("-nm",NFILE
,fnm
), 1, &nmic
, &micelle
, &micname
);
308 calc_h2order(ftp2fn(efTRX
,NFILE
,fnm
), index
, ngx
, &slDipole
, &slOrder
,
309 &slWidth
, &nslices
, top
, ePBC
, axis
, bMicel
, micelle
, nmic
,
312 h2order_plot(slDipole
, slOrder
, opt2fn("-o",NFILE
,fnm
), nslices
,
315 do_view(oenv
,opt2fn("-o",NFILE
,fnm
),"-nxy"); /* view xvgr file */