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 * Getting the Right Output Means no Artefacts in Calculating Stuff
44 int main(int argc
,char *argv
[])
48 char *fn1
,*fn2
,title
[256];
56 { efLOG
, "-f1", NULL
, ffREAD
},
57 { efLOG
, "-f2", NULL
, ffREAD
},
58 { efOUT
, "-o", "compnl", ffWRITE
},
59 { efSTX
, "-c", NULL
, ffOPTRD
}
61 #define NFILE asize(fnm)
62 static int natoms
=648;
63 static t_pargs pa
[] = {
64 { "-nat", FALSE
, etINT
, { &natoms
}, "Number of atoms" }
67 CopyRight(stderr
,argv
[0]);
68 parse_common_args(&argc
,argv
,0,TRUE
,NFILE
,fnm
,asize(pa
),pa
,0,NULL
,0,NULL
);
70 bConf
= (opt2bSet("-c",NFILE
,fnm
));
72 get_stx_coordnum (opt2fn("-c",NFILE
,fnm
),&natoms
);
73 init_t_atoms(&atoms
,natoms
,FALSE
);
75 read_stx_conf(opt2fn("-c",NFILE
,fnm
),title
,&atoms
,x
,NULL
,box
);
82 for(i
=0; (i
<natoms
); i
++) {
83 snew(mat
[0][i
],natoms
);
84 snew(mat
[1][i
],natoms
);
86 out
= ffopen(ftp2fn(efOUT
,NFILE
,fnm
),"w");
87 fn1
= opt2fn("-f1",NFILE
,fnm
);
88 fn2
= opt2fn("-f2",NFILE
,fnm
);
90 for(i
=0; (i
<2); i
++) {
91 in
= ffopen(fnm
[i
].fn
,"r");
92 fprintf(stderr
,"Reading %s\n",fnm
[i
].fn
);
93 fprintf(out
, "Reading %s\n",fnm
[i
].fn
);
94 read_nblist(in
,out
,mat
[i
],natoms
);
98 fprintf(stderr
,"Comparing Interaction Matrices\n");
99 fprintf(out
, "Comparing Interaction Matrices\n");
101 for(i
=0; (i
<natoms
); i
+=mod
) {
102 for(j
=0; (j
<natoms
); j
+=mod
) {
103 if (mat
[0][i
][j
] != mat
[1][i
][j
]) {
104 fprintf(out
,"i: %5d, j: %5d, shift[%s]: %3d, shift[%s]: %3d",
105 i
,j
,fn1
,mat
[0][i
][j
]-1,fn2
,mat
[1][i
][j
]-1);
107 pbc_dx(x
[i
],x
[j
],dx
);
108 fprintf(out
," dist: %8.3f\n",norm(dx
));
116 fprintf(out
,"There were %d mismatches\n",nmiss
);
117 fprintf(out
,"Done.\n");
119 fprintf(stderr
,"There were %d mismatches\n",nmiss
);
120 fprintf(stderr
,"Finished\n");