3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
41 #include "gmx_fatal.h"
44 void replace_atom(t_topology
*top
,int inr
,char *anm
,char *resnm
,
45 real q
,real m
,int type
)
49 atoms
= &(top
->atoms
);
51 /* Replace important properties of an atom by other properties */
52 if ((inr
< 0) || (inr
> atoms
->nr
))
53 gmx_fatal(FARGS
,"Replace_atom: inr (%d) not in %d .. %d",inr
,0,atoms
->nr
);
55 fprintf(debug
,"Replacing atom %d ... ",inr
);
56 /* Charge, mass and type */
57 atoms
->atom
[inr
].q
= atoms
->atom
[inr
].qB
= q
;
58 atoms
->atom
[inr
].m
= atoms
->atom
[inr
].mB
= m
;
59 atoms
->atom
[inr
].type
= atoms
->atom
[inr
].typeB
= type
;
62 atoms
->resinfo
[atoms
->atom
[inr
].resind
].name
= put_symtab(&top
->symtab
,resnm
);
64 atoms
->atomname
[inr
] = put_symtab(&top
->symtab
,anm
);
66 fprintf(debug
," done\n");
69 static void delete_from_interactions(t_idef
*idef
,int inr
)
75 /* Delete interactions including atom inr from lists */
76 for(i
=0; (i
<F_NRE
); i
++) {
77 nra
= interaction_function
[i
].nratoms
;
79 snew(niatoms
,idef
->il
[i
].nr
);
80 for(j
=0; (j
<idef
->il
[i
].nr
); j
+=nra
+1) {
82 for(k
=0; (k
<nra
); k
++)
83 if (idef
->il
[i
].iatoms
[j
+k
+1] == inr
)
86 /* If this does not need to be deleted, then copy it to temp array */
87 for(k
=0; (k
<nra
+1); k
++)
88 niatoms
[nnr
+k
] = idef
->il
[i
].iatoms
[j
+k
];
92 /* Copy temp array back */
93 for(j
=0; (j
<nnr
); j
++)
94 idef
->il
[i
].iatoms
[j
] = niatoms
[j
];
100 static void delete_from_block(t_block
*block
,int inr
)
102 /* Update block data structure */
105 for(i
=0; (i
<block
->nr
); i
++) {
106 for(j
=block
->index
[i
]; (j
<block
->index
[i
+1]); j
++) {
108 /* This atom has to go */
109 /* Change indices too */
110 for(i1
=i
+1; (i1
<=block
->nr
); i1
++)
117 static void delete_from_blocka(t_blocka
*block
,int inr
)
119 /* Update block data structure */
122 for(i
=0; (i
<block
->nr
); i
++) {
123 for(j
=block
->index
[i
]; (j
<block
->index
[i
+1]); j
++) {
126 /* This atom has to go */
127 for(j1
=j
; (j1
<block
->nra
-1); j1
++)
128 block
->a
[j1
] = block
->a
[j1
+1];
130 /* Change indices too */
131 for(i1
=i
+1; (i1
<=block
->nr
); i1
++)
138 static void delete_from_atoms(t_atoms
*atoms
,int inr
)
142 /* Shift the atomnames down */
143 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
144 atoms
->atomname
[i
] = atoms
->atomname
[i
+1];
146 /* Shift the atom struct down */
147 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
148 atoms
->atom
[i
] = atoms
->atom
[i
+1];
150 if (atoms
->pdbinfo
) {
151 /* Shift the pdbatom struct down */
152 for(i
=inr
; (i
<atoms
->nr
-1); i
++)
153 atoms
->pdbinfo
[i
] = atoms
->pdbinfo
[i
+1];
158 void delete_atom(t_topology
*top
,int inr
)
162 if ((inr
< 0) || (inr
>= top
->atoms
.nr
))
163 gmx_fatal(FARGS
,"Delete_atom: inr (%d) not in %d .. %d",inr
,0,
166 fprintf(debug
,"Deleting atom %d ...",inr
);
168 /* First remove bonds etc. */
169 delete_from_interactions(&top
->idef
,inr
);
170 /* Now charge groups etc. */
171 delete_from_block(&(top
->cgs
),inr
);
172 delete_from_block(&(top
->mols
),inr
);
173 delete_from_blocka(&(top
->excls
),inr
);
174 /* Now from the atoms struct */
175 delete_from_atoms(&top
->atoms
,inr
);
177 fprintf(debug
," done\n");