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