Split lines with many copyright years
[gromacs.git] / src / gromacs / topology / exclusionblocks.cpp
blobc44f53125e5f29854ec9db91db2ea54d48d1e659
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "exclusionblocks.h"
42 #include <algorithm>
44 #include "gromacs/topology/block.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/listoflists.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/utility/stringutil.h"
50 namespace gmx
53 namespace
56 //! Converts ListOfLists to a list of ExclusionBlocks
57 void listOfListsToExclusionBlocks(const ListOfLists<int>& b, gmx::ArrayRef<ExclusionBlock> b2)
59 for (gmx::index i = 0; i < b.ssize(); i++)
61 for (int jAtom : b[i])
63 b2[i].atomNumber.push_back(jAtom);
68 //! Converts a list of ExclusionBlocks to ListOfLists
69 void exclusionBlocksToListOfLists(gmx::ArrayRef<const ExclusionBlock> b2, ListOfLists<int>* b)
71 b->clear();
73 for (const auto& block : b2)
75 b->pushBack(block.atomNumber);
79 } // namespace
81 void blockaToExclusionBlocks(const t_blocka* b, gmx::ArrayRef<ExclusionBlock> b2)
83 for (int i = 0; (i < b->nr); i++)
85 for (int j = b->index[i]; (j < b->index[i + 1]); j++)
87 b2[i].atomNumber.push_back(b->a[j]);
92 void exclusionBlocksToBlocka(gmx::ArrayRef<const ExclusionBlock> b2, t_blocka* b)
94 int nra = 0;
95 int i = 0;
96 for (const auto& block : b2)
98 b->index[i] = nra;
99 int j = 0;
100 for (const auto& entry : block.atomNumber)
102 b->a[nra + j] = entry;
103 j++;
105 nra += block.nra();
106 i++;
108 /* terminate list */
109 b->index[i] = nra;
112 namespace
115 //! Counts and sorts the exclusions
116 int countAndSortExclusions(gmx::ArrayRef<ExclusionBlock> b2)
118 int nra = 0;
119 for (auto& block : b2)
121 if (block.nra() > 0)
123 /* remove double entries */
124 std::sort(block.atomNumber.begin(), block.atomNumber.end());
125 for (auto atom = block.atomNumber.begin() + 1; atom != block.atomNumber.end();)
127 GMX_RELEASE_ASSERT(atom < block.atomNumber.end(),
128 "Need to stay in range of the size of the blocks");
129 auto prev = atom - 1;
130 if (*prev == *atom)
132 atom = block.atomNumber.erase(atom);
134 else
136 ++atom;
139 nra += block.nra();
143 return nra;
146 } // namespace
148 void mergeExclusions(ListOfLists<int>* excl, gmx::ArrayRef<ExclusionBlock> b2)
150 if (b2.empty())
152 return;
154 GMX_RELEASE_ASSERT(b2.ssize() == excl->ssize(),
155 "Cannot merge exclusions for "
156 "blocks that do not describe the same number "
157 "of particles");
159 /* Convert the t_blocka entries to ExclusionBlock form */
160 listOfListsToExclusionBlocks(*excl, b2);
162 countAndSortExclusions(b2);
164 exclusionBlocksToListOfLists(b2, excl);
167 } // namespace gmx