A freeze group can now be allowed to move rigidly in some dimension by using "freezed...
[gromacs/rigid-bodies.git] / src / mdlib / wnblist.c
blob937d3283a8597caa4cae443e4a50e38167c132c7
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 * VERSION 3.2.0
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
32 * And Hey:
33 * GROwing Monsters And Cloning Shrimps
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <stdio.h>
41 #include <string.h>
42 #include "string2.h"
43 #include "force.h"
44 #include "smalloc.h"
45 #include "ns.h"
46 #include "nrnb.h"
47 #include "gmx_fatal.h"
48 #include "macros.h"
49 #include "futil.h"
50 #include "names.h"
51 #include "domdec.h"
52 #include "gmxfio.h"
54 #define header "Neighborlist:"
56 static void write_nblist(FILE *out,gmx_domdec_t *dd,t_nblist *nblist,int nDNL)
58 int i,nii,ii,j,zi,zj0,zj1,aj,zj,nj;
59 int ca1[DD_MAXZONE],np[DD_MAXZONE];
60 gmx_domdec_zones_t *dd_zones;
62 if (nblist->nri > 0) {
63 fprintf(out,"il_name: %s Solvent opt: %s\n",
64 nrnb_str(nblist->il_code),
65 enlist_names[nblist->enlist]);
66 fprintf(out,"nri: %d npair: %d\n",nblist->nri,nblist->nrj);
67 if (dd) {
68 dd_zones = domdec_zones(dd);
70 for(zi=0; zi<dd_zones->n; zi++)
71 ca1[zi] = dd->cgindex[dd_zones->cg_range[zi+1]];
72 i = 0;
73 for(zi=0; zi<dd_zones->nizone; zi++) {
74 zj0 = dd_zones->izone[zi].j0;
75 zj1 = dd_zones->izone[zi].j1;
76 for(zj=zj0; zj<zj1; zj++)
77 np[zj] = 0;
78 while(i < nblist->nri && nblist->iinr[i] < ca1[zi]) {
79 for(j=nblist->jindex[i]; (j<nblist->jindex[i+1]); j++) {
80 aj = nblist->jjnr[j];
81 zj = zj0;
82 while (aj >= ca1[zj])
83 zj++;
84 np[zj]++;
86 i++;
88 fprintf(out,"DD zone %d:",zi);
89 for(zj=zj0; zj<zj1; zj++)
90 fprintf(out," %d %d",zj,np[zj]);
91 fprintf(out,"\n");
94 if (nDNL >= 2) {
95 for(i=0; i<nblist->nri; i++) {
96 nii = 1;
97 if (nDNL >= 3 && nblist->enlist != enlistATOM_ATOM)
98 nii = 3;
99 nj = nblist->jindex[i+1] - nblist->jindex[i];
100 fprintf(out,"i: %d shift: %d gid: %d nj: %d\n",
101 ddglatnr(dd,nblist->iinr[i]),
102 nblist->shift[i],nblist->gid[i],nj);
103 for(ii=0; ii<nii; ii++) {
104 for(j=nblist->jindex[i]; (j<nblist->jindex[i+1]); j++) {
105 fprintf(out," i: %5d j: %5d\n",
106 ddglatnr(dd,nblist->iinr[i]+ii),
107 ddglatnr(dd,nblist->jjnr[j]));
112 fflush(out);
116 static void set_mat(FILE *fp,int **mat,int i0,int ni,int j0,int nj,
117 gmx_bool bSymm,int shift)
119 int i,j;
121 for(i=i0; (i<i0+ni); i++) {
122 for(j=j0; (j<j0+nj); j++) {
123 if (mat[i][j] != 0)
124 fprintf(fp,"mat[%d][%d] changing from %d to %d\n",
125 i,j,mat[i][j],shift+1);
126 mat[i][j] = shift+1;
127 if (bSymm)
128 mat[j][i] = 27-shift;
133 int read_nblist(FILE *in,FILE *fp,int **mat,int natoms,gmx_bool bSymm)
135 gmx_bool bNL;
136 char buf[256],b1[32],b2[32],solv[256],il_code[256];
137 int i,ii,j,nnbl,full,icmp,nri,isolv;
138 int iatom,nrj,nj,shift,gid,nargs,njtot=0;
140 do {
141 if (fgets2(buf,255,in) == NULL)
142 gmx_fatal(FARGS,"EOF when looking for '%s' in logfile",header);
143 } while (strstr(buf,header) == NULL);
145 do {
146 do {
147 if (fgets2(buf,255,in) == NULL)
148 return njtot;
149 } while (strstr(buf,"nri:") == NULL);
151 if (0) {
152 if ((nargs = sscanf(buf,"%*s%s%*s%s",il_code,solv)) != 2) {
153 fprintf(stderr,"Can not find the right il_code\n");
154 return njtot;
156 for(isolv=0; (isolv<esolNR); isolv++)
157 if (strstr(esol_names[isolv],solv) != NULL)
158 break;
160 if (isolv == esolNR) {
161 fprintf(stderr,"Can not read il_code or solv (nargs=%d)\n",nargs);
162 return njtot;
165 else
166 isolv = enlistATOM_ATOM;
168 /* gmx_fatal(FARGS,"Can not read il_code or solv (nargs=%d)",nargs);*/
169 if ((nargs = sscanf(buf,"%*s%d%*s%d",&nri,&nrj)) != 2)
170 gmx_fatal(FARGS,"Can not read nri or nrj (nargs=%d)",nargs);
171 for(ii=0; (ii<nri); ii++) {
172 if ((nargs = fscanf(in,"%*s%d%*s%d%*s%d%*s%d",
173 &iatom,&shift,&gid,&nj)) != 4)
174 gmx_fatal(FARGS,"Can not read iatom, shift gid or nj (nargs=%d)",nargs);
175 /* Number shifts from 1 to 27 iso 0 to 26 to distinguish uninitialized
176 * matrix elements.
178 range_check(iatom,0,natoms);
179 for(i=0; (i<nj); i++) {
180 if ((nargs = fscanf(in,"%*s%d",&j)) != 1)
181 gmx_fatal(FARGS,"Can not read j");
182 range_check(j,0,natoms);
183 switch (isolv) {
184 case enlistATOM_ATOM:
185 set_mat(fp,mat,iatom,1,j,1,bSymm,shift);
186 njtot++;
187 break;
188 case enlistSPC_ATOM:
189 set_mat(fp,mat,iatom,3,j,1,bSymm,shift);
190 njtot+=3;
191 break;
192 case enlistSPC_SPC:
193 set_mat(fp,mat,iatom,3,j,3,bSymm,shift);
194 njtot+=9;
195 break;
196 case enlistTIP4P_ATOM:
197 set_mat(fp,mat,iatom,4,j,1,bSymm,shift);
198 njtot+=4;
199 break;
200 case enlistTIP4P_TIP4P:
201 set_mat(fp,mat,iatom,4,j,4,bSymm,shift);
202 njtot+=16;
203 break;
204 gmx_incons("non-existing solvent type");
208 fprintf(fp,"nri = %d nrj = %d\n",nri,nrj);
209 } while (TRUE);
210 return -1;
213 void dump_nblist(FILE *out,t_commrec *cr,t_forcerec *fr,int nDNL)
215 #if 0
216 static FILE *fp=NULL;
217 char buf[STRLEN];
218 int n,i;
220 if (fp == NULL) {
221 if (PAR(cr)) {
222 sprintf(buf,"nlist_n%d.txt",cr->nodeid);
223 } else {
224 sprintf(buf,"nlist.txt");
226 fp = gmx_fio_fopen(buf,"w");
228 fprintf(fp,"%s\n",header);
230 for(n=0; (n<fr->nnblists); n++)
231 for(i=0; (i<eNL_NR); i++)
232 write_nblist(fp,cr->dd,&fr->nblists[n].nlist_sr[i],nDNL);
233 #endif
234 char buf[STRLEN];
235 int n,i;
237 fprintf(out,"%s\n",header);
239 for(n=0; (n<fr->nnblists); n++)
240 for(i=0; (i<eNL_NR); i++)
241 write_nblist(out,cr->dd,&fr->nblists[n].nlist_sr[i],nDNL);