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