Moved timing source to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / edittop.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/typedefs.h"
40 #include "gromacs/topology/symtab.h"
41 #include "gromacs/utility/fatalerror.h"
42 #include "gromacs/utility/smalloc.h"
43
44 void replace_atom(t_topology *top, int inr, char *anm, char *resnm,
45                   real q, real m, int type)
46 {
47     t_atoms *atoms;
48
49     atoms = &(top->atoms);
50
51     /* Replace important properties of an atom by other properties */
52     if ((inr < 0) || (inr > atoms->nr))
53     {
54         gmx_fatal(FARGS, "Replace_atom: inr (%d) not in %d .. %d", inr, 0, atoms->nr);
55     }
56     if (debug)
57     {
58         fprintf(debug, "Replacing atom %d ... ", inr);
59     }
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;
64
65     /* Residue name */
66     atoms->resinfo[atoms->atom[inr].resind].name = put_symtab(&top->symtab, resnm);
67     /* Atom name */
68     atoms->atomname[inr] = put_symtab(&top->symtab, anm);
69     if (debug)
70     {
71         fprintf(debug, " done\n");
72     }
73 }
74
75 static void delete_from_interactions(t_idef *idef, int inr)
76 {
77     int      i, j, k, nra, nnr;
78     t_iatom *niatoms;
79     gmx_bool bDel;
80
81     /* Delete interactions including atom inr from lists */
82     for (i = 0; (i < F_NRE); i++)
83     {
84         nra = interaction_function[i].nratoms;
85         nnr = 0;
86         snew(niatoms, idef->il[i].nr);
87         for (j = 0; (j < idef->il[i].nr); j += nra+1)
88         {
89             bDel = FALSE;
90             for (k = 0; (k < nra); k++)
91             {
92                 if (idef->il[i].iatoms[j+k+1] == inr)
93                 {
94                     bDel = TRUE;
95                 }
96             }
97             if (!bDel)
98             {
99                 /* If this does not need to be deleted, then copy it to temp array */
100                 for (k = 0; (k < nra+1); k++)
101                 {
102                     niatoms[nnr+k] = idef->il[i].iatoms[j+k];
103                 }
104                 nnr += nra+1;
105             }
106         }
107         /* Copy temp array back */
108         for (j = 0; (j < nnr); j++)
109         {
110             idef->il[i].iatoms[j] = niatoms[j];
111         }
112         idef->il[i].nr = nnr;
113         sfree(niatoms);
114     }
115 }
116
117 static void delete_from_block(t_block *block, int inr)
118 {
119     /* Update block data structure */
120     int i, i1, j1, j, k;
121
122     for (i = 0; (i < block->nr); i++)
123     {
124         for (j = block->index[i]; (j < block->index[i+1]); j++)
125         {
126             if (j == inr)
127             {
128                 /* This atom has to go */
129                 /* Change indices too */
130                 for (i1 = i+1; (i1 <= block->nr); i1++)
131                 {
132                     block->index[i1]--;
133                 }
134             }
135         }
136     }
137 }
138
139 static void delete_from_blocka(t_blocka *block, int inr)
140 {
141     /* Update block data structure */
142     int i, i1, j1, j, k;
143
144     for (i = 0; (i < block->nr); i++)
145     {
146         for (j = block->index[i]; (j < block->index[i+1]); j++)
147         {
148             k = block->a[j];
149             if (k == inr)
150             {
151                 /* This atom has to go */
152                 for (j1 = j; (j1 < block->nra-1); j1++)
153                 {
154                     block->a[j1] = block->a[j1+1];
155                 }
156                 block->nra--;
157                 /* Change indices too */
158                 for (i1 = i+1; (i1 <= block->nr); i1++)
159                 {
160                     block->index[i1]--;
161                 }
162             }
163         }
164     }
165 }
166
167 static void delete_from_atoms(t_atoms *atoms, int inr)
168 {
169     int i;
170
171     /* Shift the atomnames down */
172     for (i = inr; (i < atoms->nr-1); i++)
173     {
174         atoms->atomname[i] = atoms->atomname[i+1];
175     }
176
177     /* Shift the atom struct down */
178     for (i = inr; (i < atoms->nr-1); i++)
179     {
180         atoms->atom[i] = atoms->atom[i+1];
181     }
182
183     if (atoms->pdbinfo)
184     {
185         /* Shift the pdbatom struct down */
186         for (i = inr; (i < atoms->nr-1); i++)
187         {
188             atoms->pdbinfo[i] = atoms->pdbinfo[i+1];
189         }
190     }
191     atoms->nr--;
192 }
193
194 void delete_atom(t_topology *top, int inr)
195 {
196     int k;
197
198     if ((inr < 0) || (inr >= top->atoms.nr))
199     {
200         gmx_fatal(FARGS, "Delete_atom: inr (%d) not in %d .. %d", inr, 0,
201                   top->atoms.nr);
202     }
203     if (debug)
204     {
205         fprintf(debug, "Deleting atom %d ...", inr);
206     }
207
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);
216     if (debug)
217     {
218         fprintf(debug, " done\n");
219     }
220 }