Add coolquotes
[gromacs.git] / src / gromacs / topology / exclusionblocks.cpp
blobd78c2f8e2027409f7d7238ab864bc7ef0a97f864
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "exclusionblocks.h"
41 #include <algorithm>
43 #include "gromacs/topology/block.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/stringutil.h"
48 namespace gmx
51 void initExclusionBlocks(ExclusionBlocks *b2, int natom)
53 int i;
55 b2->nr = natom;
56 snew(b2->nra, b2->nr);
57 snew(b2->a, b2->nr);
58 for (i = 0; (i < b2->nr); i++)
60 b2->a[i] = nullptr;
64 void doneExclusionBlocks(ExclusionBlocks *b2)
66 int i;
68 if (b2->nr)
70 for (i = 0; (i < b2->nr); i++)
72 sfree(b2->a[i]);
74 sfree(b2->a);
75 sfree(b2->nra);
76 b2->nr = 0;
80 void blockaToExclusionBlocks(t_blocka *b, ExclusionBlocks *b2)
82 int i;
83 int j, a;
85 for (i = 0; (i < b->nr); i++)
87 for (j = b->index[i]; (j < b->index[i+1]); j++)
89 a = b->a[j];
90 srenew(b2->a[i], ++b2->nra[i]);
91 b2->a[i][b2->nra[i]-1] = a;
96 void exclusionBlocksToBlocka(ExclusionBlocks *b2, t_blocka *b)
98 int i, nra;
99 int j;
101 nra = 0;
102 for (i = 0; (i < b2->nr); i++)
104 b->index[i] = nra;
105 for (j = 0; (j < b2->nra[i]); j++)
107 b->a[nra+j] = b2->a[i][j];
109 nra += b2->nra[i];
111 /* terminate list */
112 b->index[i] = nra;
115 void mergeExclusions(t_blocka *excl, ExclusionBlocks *b2)
117 int i, k;
118 int j;
119 int nra;
121 if (!b2->nr)
123 return;
125 GMX_RELEASE_ASSERT(b2->nr == excl->nr, "Cannot merge exclusions for "
126 "blocks that do not describe the same number "
127 "of particles");
129 /* First copy all entries from excl to b2 */
130 blockaToExclusionBlocks(excl, b2);
132 /* Count and sort the exclusions */
133 nra = 0;
134 for (i = 0; (i < b2->nr); i++)
136 if (b2->nra[i] > 0)
138 /* remove double entries */
139 std::sort(b2->a[i], b2->a[i]+b2->nra[i]);
140 k = 1;
141 for (j = 1; (j < b2->nra[i]); j++)
143 if (b2->a[i][j] != b2->a[i][k-1])
145 b2->a[i][k] = b2->a[i][j];
146 k++;
149 b2->nra[i] = k;
150 nra += b2->nra[i];
153 excl->nra = nra;
154 srenew(excl->a, excl->nra);
156 exclusionBlocksToBlocka(b2, excl);
159 } // namespace gmx