Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / topology / atoms.cpp
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 "atoms.h"
40
41 #include <cstring>
42
43 #include "gromacs/topology/symtab.h"
44 #include "gromacs/utility/smalloc.h"
45
46 void init_atom(t_atoms *at)
47 {
48     at->nr        = 0;
49     at->nres      = 0;
50     at->atom      = NULL;
51     at->resinfo   = NULL;
52     at->atomname  = NULL;
53     at->atomtype  = NULL;
54     at->atomtypeB = NULL;
55     at->pdbinfo   = NULL;
56 }
57
58 void init_atomtypes(t_atomtypes *at)
59 {
60     at->nr         = 0;
61     at->radius     = NULL;
62     at->vol        = NULL;
63     at->atomnumber = NULL;
64     at->gb_radius  = NULL;
65     at->S_hct      = NULL;
66 }
67
68 void done_atom(t_atoms *at)
69 {
70     at->nr       = 0;
71     at->nres     = 0;
72     sfree(at->atom);
73     sfree(at->resinfo);
74     sfree(at->atomname);
75     sfree(at->atomtype);
76     sfree(at->atomtypeB);
77     if (at->pdbinfo)
78     {
79         sfree(at->pdbinfo);
80     }
81 }
82
83 void done_atomtypes(t_atomtypes *atype)
84 {
85     atype->nr = 0;
86     sfree(atype->radius);
87     sfree(atype->vol);
88     sfree(atype->surftens);
89     sfree(atype->atomnumber);
90     sfree(atype->gb_radius);
91     sfree(atype->S_hct);
92 }
93
94 void add_t_atoms(t_atoms *atoms, int natom_extra, int nres_extra)
95 {
96     int i;
97
98     if (natom_extra > 0)
99     {
100         srenew(atoms->atomname, atoms->nr+natom_extra);
101         srenew(atoms->atom, atoms->nr+natom_extra);
102         if (NULL != atoms->pdbinfo)
103         {
104             srenew(atoms->pdbinfo, atoms->nr+natom_extra);
105         }
106         if (NULL != atoms->atomtype)
107         {
108             srenew(atoms->atomtype, atoms->nr+natom_extra);
109         }
110         if (NULL != atoms->atomtypeB)
111         {
112             srenew(atoms->atomtypeB, atoms->nr+natom_extra);
113         }
114         for (i = atoms->nr; (i < atoms->nr+natom_extra); i++)
115         {
116             atoms->atomname[i] = NULL;
117             memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
118             if (NULL != atoms->pdbinfo)
119             {
120                 std::memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
121             }
122             if (NULL != atoms->atomtype)
123             {
124                 atoms->atomtype[i] = NULL;
125             }
126             if (NULL != atoms->atomtypeB)
127             {
128                 atoms->atomtypeB[i] = NULL;
129             }
130         }
131         atoms->nr += natom_extra;
132     }
133     if (nres_extra > 0)
134     {
135         srenew(atoms->resinfo, atoms->nres+nres_extra);
136         for (i = atoms->nres; (i < atoms->nres+nres_extra); i++)
137         {
138             std::memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
139         }
140         atoms->nres += nres_extra;
141     }
142 }
143
144 void init_t_atoms(t_atoms *atoms, int natoms, gmx_bool bPdbinfo)
145 {
146     atoms->nr   = natoms;
147     atoms->nres = 0;
148     snew(atoms->atomname, natoms);
149     atoms->atomtype  = NULL;
150     atoms->atomtypeB = NULL;
151     snew(atoms->resinfo, natoms);
152     snew(atoms->atom, natoms);
153     if (bPdbinfo)
154     {
155         snew(atoms->pdbinfo, natoms);
156     }
157     else
158     {
159         atoms->pdbinfo = NULL;
160     }
161 }
162
163 t_atoms *copy_t_atoms(t_atoms *src)
164 {
165     t_atoms *dst;
166     int      i;
167
168     snew(dst, 1);
169     init_t_atoms(dst, src->nr, (NULL != src->pdbinfo));
170     dst->nr = src->nr;
171     if (NULL != src->atomname)
172     {
173         snew(dst->atomname, src->nr);
174     }
175     if (NULL != src->atomtype)
176     {
177         snew(dst->atomtype, src->nr);
178     }
179     if (NULL != src->atomtypeB)
180     {
181         snew(dst->atomtypeB, src->nr);
182     }
183     for (i = 0; (i < src->nr); i++)
184     {
185         dst->atom[i] = src->atom[i];
186         if (NULL != src->pdbinfo)
187         {
188             dst->pdbinfo[i] = src->pdbinfo[i];
189         }
190         if (NULL != src->atomname)
191         {
192             dst->atomname[i]  = src->atomname[i];
193         }
194         if (NULL != src->atomtype)
195         {
196             dst->atomtype[i] = src->atomtype[i];
197         }
198         if (NULL != src->atomtypeB)
199         {
200             dst->atomtypeB[i] = src->atomtypeB[i];
201         }
202     }
203     dst->nres = src->nres;
204     for (i = 0; (i < src->nres); i++)
205     {
206         dst->resinfo[i] = src->resinfo[i];
207     }
208     return dst;
209 }
210
211 void t_atoms_set_resinfo(t_atoms *atoms, int atom_ind, t_symtab *symtab,
212                          const char *resname, int resnr, unsigned char ic,
213                          int chainnum, char chainid)
214 {
215     t_resinfo *ri;
216
217     ri           = &atoms->resinfo[atoms->atom[atom_ind].resind];
218     ri->name     = put_symtab(symtab, resname);
219     ri->rtp      = NULL;
220     ri->nr       = resnr;
221     ri->ic       = ic;
222     ri->chainnum = chainnum;
223     ri->chainid  = chainid;
224 }
225
226 void free_t_atoms(t_atoms *atoms, gmx_bool bFreeNames)
227 {
228     int i;
229
230     if (bFreeNames && atoms->atomname != NULL)
231     {
232         for (i = 0; i < atoms->nr; i++)
233         {
234             if (atoms->atomname[i] != NULL)
235             {
236                 sfree(*atoms->atomname[i]);
237                 *atoms->atomname[i] = NULL;
238             }
239         }
240     }
241     if (bFreeNames && atoms->resinfo != NULL)
242     {
243         for (i = 0; i < atoms->nres; i++)
244         {
245             if (atoms->resinfo[i].name != NULL)
246             {
247                 sfree(*atoms->resinfo[i].name);
248                 *atoms->resinfo[i].name = NULL;
249             }
250         }
251     }
252     if (bFreeNames && atoms->atomtype != NULL)
253     {
254         for (i = 0; i < atoms->nr; i++)
255         {
256             if (atoms->atomtype[i] != NULL)
257             {
258                 sfree(*atoms->atomtype[i]);
259                 *atoms->atomtype[i] = NULL;
260             }
261         }
262     }
263     if (bFreeNames && atoms->atomtypeB != NULL)
264     {
265         for (i = 0; i < atoms->nr; i++)
266         {
267             if (atoms->atomtypeB[i] != NULL)
268             {
269                 sfree(*atoms->atomtypeB[i]);
270                 *atoms->atomtypeB[i] = NULL;
271             }
272         }
273     }
274     sfree(atoms->atomname);
275     sfree(atoms->atomtype);
276     sfree(atoms->atomtypeB);
277     sfree(atoms->resinfo);
278     sfree(atoms->atom);
279     sfree(atoms->pdbinfo);
280     atoms->nr        = 0;
281     atoms->nres      = 0;
282     atoms->atomname  = NULL;
283     atoms->atomtype  = NULL;
284     atoms->atomtypeB = NULL;
285     atoms->resinfo   = NULL;
286     atoms->atom      = NULL;
287     atoms->pdbinfo   = NULL;
288 }