2be1f3d7533c49ed384880923ce45049decc8680
[alexxy/gromacs.git] / src / tools / edittop.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "smalloc.h"
40 #include "string2.h"
41 #include "gmx_fatal.h"
42 #include "symtab.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     gmx_fatal(FARGS,"Replace_atom: inr (%d) not in %d .. %d",inr,0,atoms->nr);
54   if (debug)
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;
60   
61   /* Residue name */
62   atoms->resinfo[atoms->atom[inr].resind].name = put_symtab(&top->symtab,resnm);
63   /* Atom name */
64   atoms->atomname[inr] = put_symtab(&top->symtab,anm);
65   if (debug)
66     fprintf(debug," done\n");
67 }
68
69 static void delete_from_interactions(t_idef *idef,int inr)
70 {
71   int  i,j,k,nra,nnr;
72   t_iatom *niatoms;
73   gmx_bool bDel;
74   
75   /* Delete interactions including atom inr from lists */
76   for(i=0; (i<F_NRE); i++) {
77     nra = interaction_function[i].nratoms;
78     nnr = 0;
79     snew(niatoms,idef->il[i].nr);
80     for(j=0; (j<idef->il[i].nr); j+=nra+1) {
81       bDel = FALSE;
82       for(k=0; (k<nra); k++)
83         if (idef->il[i].iatoms[j+k+1] == inr)
84           bDel = TRUE;
85       if (!bDel) {
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];
89         nnr+=nra+1;
90       }
91     }
92     /* Copy temp array back */
93     for(j=0; (j<nnr); j++)
94       idef->il[i].iatoms[j] = niatoms[j];
95     idef->il[i].nr = nnr;
96     sfree(niatoms);
97   }
98 }
99
100 static void delete_from_block(t_block *block,int inr) 
101 {
102   /* Update block data structure */
103   int i,i1,j1,j,k;
104   
105   for(i=0; (i<block->nr); i++) {
106     for(j=block->index[i]; (j<block->index[i+1]); j++) {
107       if (j == inr) {
108         /* This atom has to go */
109         /* Change indices too */
110         for(i1=i+1; (i1<=block->nr); i1++)
111           block->index[i1]--;
112       }
113     }
114   }
115 }
116
117 static void delete_from_blocka(t_blocka *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     for(j=block->index[i]; (j<block->index[i+1]); j++) {
124       k = block->a[j];
125       if (k == inr) {
126         /* This atom has to go */
127         for(j1=j; (j1<block->nra-1); j1++)
128           block->a[j1] = block->a[j1+1];
129         block->nra--;
130         /* Change indices too */
131         for(i1=i+1; (i1<=block->nr); i1++)
132           block->index[i1]--;
133       }
134     }
135   }
136 }
137
138 static void delete_from_atoms(t_atoms *atoms,int inr)
139 {
140   int i;
141   
142   /* Shift the atomnames down */
143   for(i=inr; (i<atoms->nr-1); i++)
144     atoms->atomname[i] = atoms->atomname[i+1];
145   
146   /* Shift the atom struct down */
147   for(i=inr; (i<atoms->nr-1); i++)
148     atoms->atom[i] = atoms->atom[i+1];
149     
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];
154   }
155   atoms->nr--;
156 }
157
158 void delete_atom(t_topology *top,int inr)
159 {
160   int k;
161   
162   if ((inr < 0) || (inr >= top->atoms.nr))
163     gmx_fatal(FARGS,"Delete_atom: inr (%d) not in %d .. %d",inr,0,
164                 top->atoms.nr);
165   if (debug)
166     fprintf(debug,"Deleting atom %d ...",inr);
167
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);
176   if (debug)
177     fprintf(debug," done\n");
178 }