9acf7ebe29680afe1e89ff8ba2887ea99228ab38
[alexxy/gromacs.git] / src / gromacs / topology / block.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 <algorithm>
40 #include "gromacs/topology/block.h"
41
42 #include "gromacs/utility/smalloc.h"
43
44 void init_block(t_block *block)
45 {
46     block->nr           = 0;
47     block->nalloc_index = 1;
48     snew(block->index, block->nalloc_index);
49     block->index[0]     = 0;
50 }
51
52 void init_blocka(t_blocka *block)
53 {
54     block->nr           = 0;
55     block->nra          = 0;
56     block->nalloc_index = 1;
57     snew(block->index, block->nalloc_index);
58     block->index[0]     = 0;
59     block->nalloc_a     = 0;
60     block->a            = NULL;
61 }
62
63 t_blocka *new_blocka(void)
64 {
65     t_blocka *block;
66
67     snew(block, 1);
68     snew(block->index, 1);
69
70     return block;
71 }
72
73 void done_block(t_block *block)
74 {
75     block->nr    = 0;
76     sfree(block->index);
77     block->nalloc_index = 0;
78 }
79
80 void done_blocka(t_blocka *block)
81 {
82     block->nr    = 0;
83     block->nra   = 0;
84     sfree(block->index);
85     sfree(block->a);
86     block->index        = NULL;
87     block->a            = NULL;
88     block->nalloc_index = 0;
89     block->nalloc_a     = 0;
90 }
91
92 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup)
93 {
94     if (bOneIndexGroup)
95     {
96         grp->nalloc_index = 2;
97         snew(grp->index, grp->nalloc_index);
98         grp->index[0] = 0;
99         grp->index[1] = natom;
100         grp->nr       = 1;
101     }
102     else
103     {
104         grp->nalloc_index = natom+1;
105         snew(grp->index, grp->nalloc_index);
106         snew(grp->index, natom+1);
107         for (int i = 0; i <= natom; ++i)
108         {
109             grp->index[i] = i;
110         }
111         grp->nr = natom;
112     }
113 }
114
115 void stupid_fill_blocka(t_blocka *grp, int natom)
116 {
117     grp->nalloc_a = natom;
118     snew(grp->a, grp->nalloc_a);
119     for (int i = 0; i < natom; ++i)
120     {
121         grp->a[i] = i;
122     }
123     grp->nra = natom;
124
125     grp->nalloc_index = natom + 1;
126     snew(grp->index, grp->nalloc_index);
127     for (int i = 0; i <= natom; ++i)
128     {
129         grp->index[i] = i;
130     }
131     grp->nr = natom;
132 }
133
134 void copy_blocka(const t_blocka *src, t_blocka *dest)
135 {
136     dest->nr           = src->nr;
137     /* Workaround for inconsistent handling of nalloc_index in
138      * other parts of the code. Often nalloc_index and nalloc_a
139      * are not set.
140      */
141     dest->nalloc_index = std::max(src->nalloc_index, dest->nr + 1);
142     snew(dest->index, dest->nalloc_index);
143     for (int i = 0; i < dest->nr+1; ++i)
144     {
145         dest->index[i] = src->index[i];
146     }
147     dest->nra      = src->nra;
148     /* See above. */
149     dest->nalloc_a = std::max(src->nalloc_a, dest->nra);
150     snew(dest->a, dest->nalloc_a);
151     for (int i = 0; i < dest->nra; ++i)
152     {
153         dest->a[i] = src->a[i];
154     }
155 }