4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * GRowing Old MAkes el Chrono Sweat
29 static char *SRCID_edittop_c
= "$Id$";
37 void replace_atom(t_topology
*top
,int inr
,char *anm
,char *resnm
,
38 real q
,real m
,int type
)
42 atoms
= &(top
->atoms
);
44 /* Replace important properties of an atom by other properties */
45 if ((inr
< 0) || (inr
> atoms
->nr
))
46 fatal_error(0,"Replace_atom: inr (%d) not in %d .. %d",inr
,0,atoms
->nr
);
48 fprintf(debug
,"Replacing atom %d ... ",inr
);
49 /* Charge, mass and type */
50 atoms
->atom
[inr
].q
= atoms
->atom
[inr
].qB
= q
;
51 atoms
->atom
[inr
].m
= atoms
->atom
[inr
].mB
= m
;
52 atoms
->atom
[inr
].type
= atoms
->atom
[inr
].typeB
= type
;
55 atoms
->resname
[atoms
->atom
[inr
].resnr
] = put_symtab(&top
->symtab
,resnm
);
57 atoms
->atomname
[inr
] = put_symtab(&top
->symtab
,anm
);
59 fprintf(debug
," done\n");
62 static void delete_from_interactions(t_idef
*idef
,int inr
)
68 /* Delete interactions including atom inr from lists */
69 for(i
=0; (i
<F_NRE
); i
++) {
70 nra
= interaction_function
[i
].nratoms
;
72 snew(niatoms
,idef
->il
[i
].nr
);
73 for(j
=0; (j
<idef
->il
[i
].nr
); j
+=nra
+1) {
75 for(k
=0; (k
<nra
); k
++)
76 if (idef
->il
[i
].iatoms
[j
+k
+1] == inr
)
79 /* If this does not need to be deleted, then copy it to temp array */
80 for(k
=0; (k
<nra
+1); k
++)
81 niatoms
[nnr
+k
] = idef
->il
[i
].iatoms
[j
+k
];
85 /* Copy temp array back */
86 for(j
=0; (j
<nnr
); j
++)
87 idef
->il
[i
].iatoms
[j
] = niatoms
[j
];
91 /* Reduce multinr if necessary */
92 for(j
=0; (j
<MAXPROC
); j
++)
93 if (idef
->il
[i
].multinr
[j
] >= nnr
)
94 idef
->il
[i
].multinr
[j
] = nnr
;
98 static void delete_from_block(t_block
*block
,int inr
)
100 /* Update block data structure */
103 for(i
=0; (i
<block
->nr
); i
++) {
104 for(j
=block
->index
[i
]; (j
<block
->index
[i
+1]); j
++) {
107 /* This atom has to go */
108 for(j1
=j
; (j1
<block
->nra
-1); j1
++)
109 block
->a
[j1
] = block
->a
[j1
+1];
111 /* Change indices too */
112 for(i1
=i
+1; (i1
<=block
->nr
); i1
++)
119 static void delete_from_atoms(t_atoms
*atoms
,int inr
)
123 /* Delete inr as an exclusion from other atoms */
124 delete_from_block(&(atoms
->excl
),inr
);
125 /* Now delete the exclusions with inr as i atom */
126 ind0
= atoms
->excl
.index
[inr
];
127 nrei
= atoms
->excl
.index
[inr
+1]-ind0
;
128 for(i
=ind0
+nrei
; (i
<atoms
->excl
.nra
); i
++)
129 atoms
->excl
.a
[i
-nrei
] = atoms
->excl
.a
[i
];
130 atoms
->excl
.nra
-= nrei
;
131 for(i
=inr
; (i
<atoms
->excl
.nr
); i
++)
132 atoms
->excl
.index
[i
] = atoms
->excl
.index
[i
+1] - nrei
;
134 assert(atoms
->excl
.index
[atoms
->excl
.nr
] == atoms
->excl
.nra
);
136 /* Shift the atomnames down */
137 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
138 atoms
->atomname
[i
] = atoms
->atomname
[i
+1];
140 /* Shift the atom struct down */
141 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
142 atoms
->atom
[i
] = atoms
->atom
[i
+1];
144 if (atoms
->pdbinfo
) {
145 /* Shift the pdbatom struct down */
146 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
147 atoms
->pdbinfo
[i
] = atoms
->pdbinfo
[i
+1];
152 void delete_atom(t_topology
*top
,int inr
)
156 if ((inr
< 0) || (inr
>= top
->atoms
.nr
))
157 fatal_error(0,"Delete_atom: inr (%d) not in %d .. %d",inr
,0,
160 fprintf(debug
,"Deleting atom %d ...",inr
);
162 /* First remove bonds etc. */
163 delete_from_interactions(&top
->idef
,inr
);
164 /* Now charge groups etc. */
165 for(k
=0; (k
<ebNR
); k
++)
166 delete_from_block(&(top
->blocks
[k
]),inr
);
167 /* Now from the atoms struct */
168 delete_from_atoms(&top
->atoms
,inr
);
170 fprintf(debug
," done\n");