4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
30 * Gromacs Runs On Most of All Computer Systems
33 /* This file is completely threadsafe - keep it that way! */
44 #define header "Neighborlist:"
46 static void write_nblist(FILE *out
,t_nblist
*nblist
)
48 int i
,j
,j0
,k
,i_atom
,jid
,nj
;
49 fprintf(out
,"il_code: %d solvent: %d\n",nblist
->il_code
,nblist
->solvent
);
51 fprintf(out
,"nri: %d nrj: %d\n",nblist
->nri
,nblist
->nrj
);
52 for(i
=0; i
<nblist
->nri
; i
++) {
53 nj
= nblist
->jindex
[i
+1] - nblist
->jindex
[i
];
54 fprintf(out
,"i: %d shift: %d gid: %d nj: %d\n",
55 nblist
->iinr
[i
],nblist
->shift
[i
],nblist
->gid
[i
],nj
);
56 for(j
=nblist
->jindex
[i
]; (j
<nblist
->jindex
[i
+1]); j
++)
57 fprintf(out
," j: %d\n",nblist
->jjnr
[j
]);
62 void read_nblist(FILE *in
,FILE *log
,int **mat
,int natoms
)
65 char buf
[256],b1
[32],b2
[32];
66 int i
,ii
,j
,nnbl
,full
,icmp
,nri
,il_code
,solv
;
67 int iatom
,nrj
,nj
,shift
,gid
;
70 if (fgets2(buf
,255,in
) == NULL
)
71 fatal_error(0,"EOF when looking for '%s' in logfile",header
);
72 } while (strstr(buf
,header
) == NULL
);
75 if (fscanf(in
,"%*s%d%*s%d",&il_code
,&solv
) != 2)
77 if (fscanf(in
,"%*s%d%*s%d",&nri
,&nrj
) != 2)
78 fatal_error(0,"Not enough arguments read line %d",__LINE__
);
79 for(ii
=0; (ii
<nri
); ii
++) {
80 if (fscanf(in
,"%*s%d%*s%d%*s%d%*s%d",&iatom
,&gid
,&shift
,&nj
) != 4)
81 fatal_error(0,"Not enough arguments read line %d",__LINE__
);
82 /* Number shifts from 1 to 27 iso 0 to 26, to distinguish uninitialized
86 if ((iatom
< 0) || (iatom
>= natoms
))
87 fatal_error(0,"iatom = %d (max %d)\n",iatom
,natoms
);
89 for(i
=0; (i
<nj
); i
++) {
90 if (fscanf(in
,"%*s%d",&j
) != 1)
91 fatal_error(0,"Not enough arguments read line %d",__LINE__
);
92 if ((j
< 0) || (j
>= natoms
))
93 fatal_error(0,"iatom = %d (max %d)\n",j
,natoms
);
94 if (mat
[iatom
][j
] != 0)
95 fprintf(log
,"mat[%d][%d] changing from %d to %d\n",
96 i
,j
,mat
[iatom
][j
],shift
);
97 mat
[iatom
][j
] = shift
;
100 fprintf(log
,"nri = %d nrj = %d\n",nri
,nrj
);
104 void dump_nblist(FILE *out
,t_forcerec
*fr
,int nDNL
)
108 fprintf(out
,"%s\n",header
);
110 for(i
=0; (i
<eNL_NR
); i
++)
111 write_nblist(out
,&fr
->nlist_sr
[i
]);