Removed include of types/ifunc.h from typedefs.h
[gromacs.git] / src / gromacs / gmxana / edittop.cpp
blobd7b55b6baf6f99c4c661dbc09fbe61364860d328
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "gromacs/legacyheaders/typedefs.h"
40 #include "gromacs/legacyheaders/types/ifunc.h"
41 #include "gromacs/topology/symtab.h"
42 #include "gromacs/utility/fatalerror.h"
43 #include "gromacs/utility/smalloc.h"
45 void replace_atom(t_topology *top, int inr, char *anm, char *resnm,
46 real q, real m, int type)
48 t_atoms *atoms;
50 atoms = &(top->atoms);
52 /* Replace important properties of an atom by other properties */
53 if ((inr < 0) || (inr > atoms->nr))
55 gmx_fatal(FARGS, "Replace_atom: inr (%d) not in %d .. %d", inr, 0, atoms->nr);
57 if (debug)
59 fprintf(debug, "Replacing atom %d ... ", inr);
61 /* Charge, mass and type */
62 atoms->atom[inr].q = atoms->atom[inr].qB = q;
63 atoms->atom[inr].m = atoms->atom[inr].mB = m;
64 atoms->atom[inr].type = atoms->atom[inr].typeB = type;
66 /* Residue name */
67 atoms->resinfo[atoms->atom[inr].resind].name = put_symtab(&top->symtab, resnm);
68 /* Atom name */
69 atoms->atomname[inr] = put_symtab(&top->symtab, anm);
70 if (debug)
72 fprintf(debug, " done\n");
76 static void delete_from_interactions(t_idef *idef, int inr)
78 int i, j, k, nra, nnr;
79 t_iatom *niatoms;
80 gmx_bool bDel;
82 /* Delete interactions including atom inr from lists */
83 for (i = 0; (i < F_NRE); i++)
85 nra = interaction_function[i].nratoms;
86 nnr = 0;
87 snew(niatoms, idef->il[i].nr);
88 for (j = 0; (j < idef->il[i].nr); j += nra+1)
90 bDel = FALSE;
91 for (k = 0; (k < nra); k++)
93 if (idef->il[i].iatoms[j+k+1] == inr)
95 bDel = TRUE;
98 if (!bDel)
100 /* If this does not need to be deleted, then copy it to temp array */
101 for (k = 0; (k < nra+1); k++)
103 niatoms[nnr+k] = idef->il[i].iatoms[j+k];
105 nnr += nra+1;
108 /* Copy temp array back */
109 for (j = 0; (j < nnr); j++)
111 idef->il[i].iatoms[j] = niatoms[j];
113 idef->il[i].nr = nnr;
114 sfree(niatoms);
118 static void delete_from_block(t_block *block, int inr)
120 /* Update block data structure */
121 int i, i1, j;
123 for (i = 0; (i < block->nr); i++)
125 for (j = block->index[i]; (j < block->index[i+1]); j++)
127 if (j == inr)
129 /* This atom has to go */
130 /* Change indices too */
131 for (i1 = i+1; (i1 <= block->nr); i1++)
133 block->index[i1]--;
140 static void delete_from_blocka(t_blocka *block, int inr)
142 /* Update block data structure */
143 int i, i1, j1, j, k;
145 for (i = 0; (i < block->nr); i++)
147 for (j = block->index[i]; (j < block->index[i+1]); j++)
149 k = block->a[j];
150 if (k == inr)
152 /* This atom has to go */
153 for (j1 = j; (j1 < block->nra-1); j1++)
155 block->a[j1] = block->a[j1+1];
157 block->nra--;
158 /* Change indices too */
159 for (i1 = i+1; (i1 <= block->nr); i1++)
161 block->index[i1]--;
168 static void delete_from_atoms(t_atoms *atoms, int inr)
170 int i;
172 /* Shift the atomnames down */
173 for (i = inr; (i < atoms->nr-1); i++)
175 atoms->atomname[i] = atoms->atomname[i+1];
178 /* Shift the atom struct down */
179 for (i = inr; (i < atoms->nr-1); i++)
181 atoms->atom[i] = atoms->atom[i+1];
184 if (atoms->pdbinfo)
186 /* Shift the pdbatom struct down */
187 for (i = inr; (i < atoms->nr-1); i++)
189 atoms->pdbinfo[i] = atoms->pdbinfo[i+1];
192 atoms->nr--;
195 void delete_atom(t_topology *top, int inr)
197 if ((inr < 0) || (inr >= top->atoms.nr))
199 gmx_fatal(FARGS, "Delete_atom: inr (%d) not in %d .. %d", inr, 0,
200 top->atoms.nr);
202 if (debug)
204 fprintf(debug, "Deleting atom %d ...", inr);
207 /* First remove bonds etc. */
208 delete_from_interactions(&top->idef, inr);
209 /* Now charge groups etc. */
210 delete_from_block(&(top->cgs), inr);
211 delete_from_block(&(top->mols), inr);
212 delete_from_blocka(&(top->excls), inr);
213 /* Now from the atoms struct */
214 delete_from_atoms(&top->atoms, inr);
215 if (debug)
217 fprintf(debug, " done\n");