1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
|
/*
* $Id: stat.c,v 1.29.2.3 2007/09/20 06:35:41 spoel Exp $
*
* This source code is part of
*
* G R O M A C S
*
* GROningen MAchine for Chemical Simulations
*
* VERSION 3.3.2
* Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
* Copyright (c) 1991-2000, University of Groningen, The Netherlands.
* Copyright (c) 2001-2007, The GROMACS development team,
* check out http://www.gromacs.org for more information.
* This program is free software; you can redistribute it and/or
* modify it under the terms of the GNU General Public License
* as published by the Free Software Foundation; either version 2
* of the License, or (at your option) any later version.
*
* If you want to redistribute modifications, please consider that
* scientific software is very special. Version control is crucial -
* bugs must be traceable. We will be happy to consider code for
* inclusion in the official distribution, but derived work must not
* be called official GROMACS. Details are found in the README & COPYING
* files - if they are missing, get the official version at www.gromacs.org.
*
* To help us fund GROMACS development, we humbly ask that you cite
* the papers on the package - you can find them in the top README file.
*
* For more info, check our website at http://www.gromacs.org
*
* And Hey:
* Groningen Machine for Chemical Simulation
*/
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#include <string.h>
#include <stdio.h>
#include "typedefs.h"
#include "sysstuff.h"
#include "gmx_fatal.h"
#include "network.h"
#include "txtdump.h"
#include "names.h"
#include "physics.h"
#include "vec.h"
#include "maths.h"
#include "mvdata.h"
#include "main.h"
#include "force.h"
#include "nrnb.h"
#include "vcm.h"
#include "smalloc.h"
#include "futil.h"
#include "network.h"
#include "rbin.h"
#include "tgroup.h"
#include "xtcio.h"
#include "gmxfio.h"
#include "trnio.h"
#include "statutil.h"
#include "hdf5_simresults.h"
/* STRUCT: SimResultsBond
*
* A copy of the ne1::SimResultsBond defined in SimResultsDataStore.h
*/
typedef struct SimResultsBond {
unsigned int atomId_1, atomId_2;
float order;
} SimResultsBond;
void global_stat(FILE *log,
t_commrec *cr,real ener[],
tensor fvir,tensor svir,
t_grpopts *opts,t_groups *grps,
t_nrnb *mynrnb,t_nrnb nrnb[],
t_vcm *vcm,real *terminate)
{
static t_bin *rb=NULL;
static int *itc;
int iterminate,imu,ie,ifv,isv,idedl,icm,imass,ica;
int icj=-1,ici=-1,icx=-1;
int in[MAXNODES];
int inn[egNR];
int j;
if (rb==NULL) {
rb=mk_bin();
snew(itc,opts->ngtc);
}
else
reset_bin(rb);
/* Reset nrnb stuff */
for(j=0; (j<cr->nnodes); j++)
init_nrnb(&(nrnb[j]));
cp_nrnb(&(nrnb[cr->nodeid]),mynrnb);
/* This routine copies all the data to be summed to one big buffer
* using the t_bin struct.
*/
where();
ie = add_binr(log,rb,F_NRE,ener);
where();
ifv = add_binr(log,rb,DIM*DIM,fvir[0]);
where();
isv = add_binr(log,rb,DIM*DIM,svir[0]);
where();
for(j=0; (j<cr->nnodes); j++)
in[j] = add_bind(log,rb,eNRNB,nrnb[j].n);
where();
for(j=0; (j<opts->ngtc); j++)
itc[j]=add_binr(log,rb,DIM*DIM,grps->tcstat[j].ekinh[0]);
where();
idedl = add_binr(log,rb,1,&(grps->dekindl));
where();
for(j=0; (j<egNR); j++)
inn[j]=add_binr(log,rb,grps->estat.nn,grps->estat.ee[j]);
where();
icm = add_binr(log,rb,DIM*vcm->nr,vcm->group_p[0]);
where();
imass = add_binr(log,rb,vcm->nr,vcm->group_mass);
where();
if (vcm->mode == ecmANGULAR) {
icj = add_binr(log,rb,DIM*vcm->nr,vcm->group_j[0]);
where();
icx = add_binr(log,rb,DIM*vcm->nr,vcm->group_x[0]);
where();
ici = add_binr(log,rb,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
where();
}
ica = add_binr(log,rb,1,&(grps->cosacc.mvcos));
where();
iterminate = add_binr(log,rb,1,terminate);
/* Global sum it all */
sum_bin(rb,cr);
where();
/* Extract all the data locally */
extract_binr(rb,ie ,F_NRE,ener);
extract_binr(rb,ifv ,DIM*DIM,fvir[0]);
extract_binr(rb,isv ,DIM*DIM,svir[0]);
for(j=0; (j<cr->nnodes); j++)
extract_bind(rb,in[j],eNRNB,nrnb[j].n);
for(j=0; (j<opts->ngtc); j++)
extract_binr(rb,itc[j],DIM*DIM,grps->tcstat[j].ekinh[0]);
extract_binr(rb,idedl,1,&(grps->dekindl));
for(j=0; (j<egNR); j++)
extract_binr(rb,inn[j],grps->estat.nn,grps->estat.ee[j]);
extract_binr(rb,icm,DIM*vcm->nr,vcm->group_p[0]);
where();
extract_binr(rb,imass,vcm->nr,vcm->group_mass);
where();
if (vcm->mode == ecmANGULAR) {
extract_binr(rb,icj,DIM*vcm->nr,vcm->group_j[0]);
where();
extract_binr(rb,icx,DIM*vcm->nr,vcm->group_x[0]);
where();
extract_binr(rb,ici,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
where();
}
extract_binr(rb,ica,1,&(grps->cosacc.mvcos));
where();
extract_binr(rb,iterminate,1,terminate);
where();
/* Small hack for temp only */
ener[F_TEMP]/=cr->nnodes;
}
int do_per_step(int step,int nstep)
{
if (nstep != 0)
return ((step % nstep)==0);
else
return 0;
}
static void moveit(FILE *log,
int left,int right,char *s,rvec xx[],t_nsborder *nsb)
{
if (!xx)
return;
move_rvecs(log,FALSE,FALSE,left,right,xx,NULL,nsb->nnodes-1,nsb,NULL);
}
int write_traj(FILE *log,t_commrec *cr,
char *traj,t_nsborder *nsb,
int step,real t,real lambda,t_nrnb nrnb[],
int natoms,rvec *xx,rvec *vv,rvec *ff,matrix box,
t_topology *top)
{
static int fp=-1;
int initializeHDF5_Frame = 0;
if ((fp == -1) && MASTER(cr)) {
#ifdef DEBUG
fprintf(log,"Going to open trajectory file: %s\n",traj);
#endif
fp = open_trn(traj,"w");
if (fn2ftp(traj) == efNH5)
initializeHDF5_Frame = 1; // Write bonds to first frame.
}
#define MX(xvf) moveit(log,cr->left,cr->right,#xvf,xvf,nsb)
if (cr->nnodes > 1) {
MX(xx);
MX(vv);
MX(ff);
}
if ((xx || vv || ff) && MASTER(cr)) {
fwrite_trn(fp,step,t,lambda,box,natoms,xx,vv,ff);
// If we're writing to a Nanorex HDF5 datastore, write bonds and atom types
// to the first frame.
//
if (initializeHDF5_Frame) {
printf(">>> stat.c:write_traj: hdf5 datastore opened - writing topo\n");
// Write atoms
//
// One pass through to determine the atom count for the non-solvent
// structure. The solvent residues will always be appended to the real
// structure.
//
unsigned int atomIndex = 0;
unsigned int atomCount = top->atoms.nr;
for (atomIndex = 0; atomIndex < atomCount; atomIndex++) {
if (strcmp((*(top->atoms.resname[top->atoms.atom[atomIndex].resnr])),
"SOL") == 0) {
atomCount = atomIndex;
break;
}
}
// Now we know the true atom count
unsigned int* atomIds =
(unsigned int*)malloc(atomCount*sizeof(unsigned int));
unsigned int* atomicNumbers =
(unsigned int*)malloc(atomCount*sizeof(unsigned int));
for (atomIndex = 0; atomIndex < atomCount; atomIndex++) {
atomIds[atomIndex] = atomIndex;
// NOTE: Atomic numbers are specified in the .itp files, just before
// the mass column.
atomicNumbers[atomIndex] =
top->atomtypes.atomnumber[top->atoms.atom[atomIndex].type];
}
addHDF5atomIds(atomIds, atomCount);
free(atomIds);
addHDF5atomicNumbers(atomicNumbers, atomCount);
free(atomicNumbers);
// Write bonds. (From src/kernel/gmxcheck.c:chk_bonds().)
//
// One loop to determine the number of bonds
//
unsigned int bondCount = 0;
int ftype, k, ai, aj;
for (ftype = 0; ftype <= F_SHAKE; ftype++) {
if ((interaction_function[ftype].flags & IF_CHEMBOND) ==
IF_CHEMBOND) {
bondCount += top->idef.il[ftype].nr / 3;
}
}
// Now create the bonds array
//
unsigned int bondIndex = 0;
SimResultsBond bond;
void* bonds = (void*)malloc(bondCount*sizeof(SimResultsBond));
for (ftype = 0; ftype <= F_SHAKE; ftype++) {
switch (ftype) {
case F_BONDS:
case F_G96BONDS:
case F_MORSE:
case F_CUBICBONDS:
case F_SHAKE:
if ((interaction_function[ftype].flags & IF_CHEMBOND) ==
IF_CHEMBOND) {
for (k = 0; k < top->idef.il[ftype].nr; ) {
k++; // "type"
bond.order = 0;
bond.atomId_1 = top->idef.il[ftype].iatoms[k++];
bond.atomId_2 = top->idef.il[ftype].iatoms[k++];
if (bondIndex >= bondCount) {
printf(">>> src/gmxlib/stat.c:write_traj: miscalculated the number of bonds - bondCount=%d\n",
bondCount);
gmx_fatal(FARGS,"HDF5 error");
}
((SimResultsBond*)bonds)[bondIndex] = bond;
bondIndex++;
}
}
break;
}
}
bondCount = bondIndex;
addHDF5bonds(bonds, bondCount);
free(bonds);
}
fio_flush(fp);
}
return fp;
}
/* XDR stuff for compressed trajectories */
static int xd;
void write_xtc_traj(FILE *log,t_commrec *cr,
char *xtc_traj,t_nsborder *nsb,t_mdatoms *md,
int step,real t,rvec *xx,matrix box,real prec)
{
static bool bFirst=TRUE;
static rvec *x_sel;
static int natoms;
int i,j;
if ((bFirst) && MASTER(cr)) {
#ifdef DEBUG
fprintf(log,"Going to open compressed trajectory file: %s\n",xtc_traj);
#endif
xd = open_xtc(xtc_traj,"w");
/* Count the number of atoms in the selection */
natoms=0;
for(i=0; (i<md->nr); i++)
if (md->cXTC[i] == 0)
natoms++;
fprintf(log,"There are %d atoms in your xtc output selection\n",natoms);
if (natoms != md->nr)
snew(x_sel,natoms);
bFirst=FALSE;
}
if (cr->nnodes > 1) {
MX(xx);
}
if ((xx) && MASTER(cr)) {
if (natoms == md->nr)
x_sel = xx;
else {
/* We need to copy everything into a temp array */
for(i=j=0; (i<md->nr); i++) {
if (md->cXTC[i] == 0) {
copy_rvec(xx[i],x_sel[j]);
j++;
}
}
}
if (write_xtc(xd,natoms,step,t,box,x_sel,prec) == 0)
gmx_fatal(FARGS,"XTC error");
}
}
void close_xtc_traj(void)
{
close_xtc(xd);
}
|