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))
54 gmx_fatal(FARGS, "Replace_atom: inr (%d) not in %d .. %d", inr, 0, atoms->nr);
58 fprintf(debug, "Replacing atom %d ... ", inr);
60 /* Charge, mass and type */
61 atoms->atom[inr].q = atoms->atom[inr].qB = q;
62 atoms->atom[inr].m = atoms->atom[inr].mB = m;
63 atoms->atom[inr].type = atoms->atom[inr].typeB = type;
66 atoms->resinfo[atoms->atom[inr].resind].name = put_symtab(&top->symtab, resnm);
68 atoms->atomname[inr] = put_symtab(&top->symtab, anm);
71 fprintf(debug, " done\n");
75 static void delete_from_interactions(t_idef *idef, int inr)
77 int i, j, k, nra, nnr;
81 /* Delete interactions including atom inr from lists */
82 for (i = 0; (i < F_NRE); i++)
84 nra = interaction_function[i].nratoms;
86 snew(niatoms, idef->il[i].nr);
87 for (j = 0; (j < idef->il[i].nr); j += nra+1)
90 for (k = 0; (k < nra); k++)
92 if (idef->il[i].iatoms[j+k+1] == inr)
99 /* If this does not need to be deleted, then copy it to temp array */
100 for (k = 0; (k < nra+1); k++)
102 niatoms[nnr+k] = idef->il[i].iatoms[j+k];
107 /* Copy temp array back */
108 for (j = 0; (j < nnr); j++)
110 idef->il[i].iatoms[j] = niatoms[j];
112 idef->il[i].nr = nnr;
117 static void delete_from_block(t_block *block, int inr)
119 /* Update block data structure */
122 for (i = 0; (i < block->nr); i++)
124 for (j = block->index[i]; (j < block->index[i+1]); j++)
128 /* This atom has to go */
129 /* Change indices too */
130 for (i1 = i+1; (i1 <= block->nr); i1++)
139 static void delete_from_blocka(t_blocka *block, int inr)
141 /* Update block data structure */
144 for (i = 0; (i < block->nr); i++)
146 for (j = block->index[i]; (j < block->index[i+1]); j++)
151 /* This atom has to go */
152 for (j1 = j; (j1 < block->nra-1); j1++)
154 block->a[j1] = block->a[j1+1];
157 /* Change indices too */
158 for (i1 = i+1; (i1 <= block->nr); i1++)
167 static void delete_from_atoms(t_atoms *atoms, int inr)
171 /* Shift the atomnames down */
172 for (i = inr; (i < atoms->nr-1); i++)
174 atoms->atomname[i] = atoms->atomname[i+1];
177 /* Shift the atom struct down */
178 for (i = inr; (i < atoms->nr-1); i++)
180 atoms->atom[i] = atoms->atom[i+1];
185 /* Shift the pdbatom struct down */
186 for (i = inr; (i < atoms->nr-1); i++)
188 atoms->pdbinfo[i] = atoms->pdbinfo[i+1];
194 void delete_atom(t_topology *top, int inr)
198 if ((inr < 0) || (inr >= top->atoms.nr))
200 gmx_fatal(FARGS, "Delete_atom: inr (%d) not in %d .. %d", inr, 0,
205 fprintf(debug, "Deleting atom %d ...", inr);
208 /* First remove bonds etc. */
209 delete_from_interactions(&top->idef, inr);
210 /* Now charge groups etc. */
211 delete_from_block(&(top->cgs), inr);
212 delete_from_block(&(top->mols), inr);
213 delete_from_blocka(&(top->excls), inr);
214 /* Now from the atoms struct */
215 delete_from_atoms(&top->atoms, inr);
218 fprintf(debug, " done\n");