7c1cc84357f905775b0f3f38749d84ea473567a1
[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,2015,2017,2018,2019, 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 "block.h"
40
41 #include <cstdio>
42
43 #include <algorithm>
44
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/txtdump.h"
47
48 void gmx::RangePartitioning::setAllBlocksSizeOne(int numBlocksToSet)
49 {
50     if (!allBlocksHaveSizeOne())
51     {
52         clear();
53     }
54     if (numBlocksToSet < numBlocks())
55     {
56         index_.resize(numBlocksToSet + 1);
57     }
58     else if (numBlocksToSet > numBlocks())
59     {
60         for (int b = numBlocks(); b < numBlocksToSet; b++)
61         {
62             appendBlock(1);
63         }
64     }
65 }
66
67 void init_block(t_block* block)
68 {
69     block->nr           = 0;
70     block->nalloc_index = 1;
71     snew(block->index, block->nalloc_index);
72     block->index[0] = 0;
73 }
74
75 void init_block_null(t_block* block)
76 {
77     block->nr           = 0;
78     block->nalloc_index = 0;
79     block->index        = nullptr;
80 }
81
82 void init_blocka(t_blocka* block)
83 {
84     block->nr           = 0;
85     block->nra          = 0;
86     block->nalloc_index = 1;
87     snew(block->index, block->nalloc_index);
88     block->index[0] = 0;
89     block->nalloc_a = 0;
90     block->a        = nullptr;
91 }
92
93 void init_blocka_null(t_blocka* block)
94 {
95     block->nr           = 0;
96     block->nra          = 0;
97     block->nalloc_index = 0;
98     block->index        = nullptr;
99     block->nalloc_a     = 0;
100     block->a            = nullptr;
101 }
102
103 t_blocka* new_blocka()
104 {
105     t_blocka* block;
106
107     snew(block, 1);
108     snew(block->index, 1);
109
110     return block;
111 }
112
113 void done_block(t_block* block)
114 {
115     block->nr = 0;
116     sfree(block->index);
117     block->index        = nullptr;
118     block->nalloc_index = 0;
119 }
120
121 void done_blocka(t_blocka* block)
122 {
123     block->nr  = 0;
124     block->nra = 0;
125     sfree(block->index);
126     sfree(block->a);
127     block->index        = nullptr;
128     block->a            = nullptr;
129     block->nalloc_index = 0;
130     block->nalloc_a     = 0;
131 }
132
133 void stupid_fill_block(t_block* grp, int natom, gmx_bool bOneIndexGroup)
134 {
135     if (bOneIndexGroup)
136     {
137         grp->nalloc_index = 2;
138         srenew(grp->index, grp->nalloc_index);
139         grp->index[0] = 0;
140         grp->index[1] = natom;
141         grp->nr       = 1;
142     }
143     else
144     {
145         grp->nalloc_index = natom + 1;
146         srenew(grp->index, grp->nalloc_index);
147         for (int i = 0; i <= natom; ++i)
148         {
149             grp->index[i] = i;
150         }
151         grp->nr = natom;
152     }
153 }
154
155 void stupid_fill_blocka(t_blocka* grp, int natom)
156 {
157     grp->nalloc_a = natom;
158     snew(grp->a, grp->nalloc_a);
159     for (int i = 0; i < natom; ++i)
160     {
161         grp->a[i] = i;
162     }
163     grp->nra = natom;
164
165     grp->nalloc_index = natom + 1;
166     snew(grp->index, grp->nalloc_index);
167     for (int i = 0; i <= natom; ++i)
168     {
169         grp->index[i] = i;
170     }
171     grp->nr = natom;
172 }
173
174 void copy_blocka(const t_blocka* src, t_blocka* dest)
175 {
176     dest->nr = src->nr;
177     /* Workaround for inconsistent handling of nalloc_index in
178      * other parts of the code. Often nalloc_index and nalloc_a
179      * are not set.
180      */
181     dest->nalloc_index = std::max(src->nalloc_index, dest->nr + 1);
182     snew(dest->index, dest->nalloc_index);
183     for (int i = 0; i < dest->nr + 1; ++i)
184     {
185         dest->index[i] = src->index[i];
186     }
187     dest->nra = src->nra;
188     /* See above. */
189     dest->nalloc_a = std::max(src->nalloc_a, dest->nra);
190     snew(dest->a, dest->nalloc_a);
191     for (int i = 0; i < dest->nra; ++i)
192     {
193         dest->a[i] = src->a[i];
194     }
195 }
196
197 static int pr_block_title(FILE* fp, int indent, const char* title, const t_block* block)
198 {
199     if (available(fp, block, indent, title))
200     {
201         indent = pr_title(fp, indent, title);
202         pr_indent(fp, indent);
203         fprintf(fp, "nr=%d\n", block->nr);
204     }
205     return indent;
206 }
207
208 static int pr_blocka_title(FILE* fp, int indent, const char* title, const t_blocka* block)
209 {
210     if (available(fp, block, indent, title))
211     {
212         indent = pr_title(fp, indent, title);
213         pr_indent(fp, indent);
214         fprintf(fp, "nr=%d\n", block->nr);
215         pr_indent(fp, indent);
216         fprintf(fp, "nra=%d\n", block->nra);
217     }
218     return indent;
219 }
220
221 static void low_pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers)
222 {
223     int i;
224
225     if (available(fp, block, indent, title))
226     {
227         indent = pr_blocka_title(fp, indent, title, block);
228         for (i = 0; i <= block->nr; i++)
229         {
230             pr_indent(fp, indent + INDENT);
231             fprintf(fp, "%s->index[%d]=%d\n", title, bShowNumbers ? i : -1, block->index[i]);
232         }
233         for (i = 0; i < block->nra; i++)
234         {
235             pr_indent(fp, indent + INDENT);
236             fprintf(fp, "%s->a[%d]=%d\n", title, bShowNumbers ? i : -1, block->a[i]);
237         }
238     }
239 }
240
241 void pr_block(FILE* fp, int indent, const char* title, const t_block* block, gmx_bool bShowNumbers)
242 {
243     int i, start;
244
245     if (available(fp, block, indent, title))
246     {
247         indent = pr_block_title(fp, indent, title, block);
248         start  = 0;
249         if (block->index[start] != 0)
250         {
251             fprintf(fp, "block->index[%d] should be 0\n", start);
252         }
253         else
254         {
255             for (i = 0; i < block->nr; i++)
256             {
257                 int end = block->index[i + 1];
258                 pr_indent(fp, indent);
259                 if (end <= start)
260                 {
261                     fprintf(fp, "%s[%d]={}\n", title, i);
262                 }
263                 else
264                 {
265                     fprintf(fp, "%s[%d]={%d..%d}\n", title, bShowNumbers ? i : -1,
266                             bShowNumbers ? start : -1, bShowNumbers ? end - 1 : -1);
267                 }
268                 start = end;
269             }
270         }
271     }
272 }
273
274 void pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers)
275 {
276     int i, j, ok, size, start, end;
277
278     if (available(fp, block, indent, title))
279     {
280         indent = pr_blocka_title(fp, indent, title, block);
281         start  = 0;
282         end    = start;
283         if ((ok = static_cast<int>(block->index[start] == 0)) == 0)
284         {
285             fprintf(fp, "block->index[%d] should be 0\n", start);
286         }
287         else
288         {
289             for (i = 0; i < block->nr; i++)
290             {
291                 end  = block->index[i + 1];
292                 size = pr_indent(fp, indent);
293                 if (end <= start)
294                 {
295                     size += fprintf(fp, "%s[%d]={", title, i);
296                 }
297                 else
298                 {
299                     size += fprintf(fp, "%s[%d][%d..%d]={", title, bShowNumbers ? i : -1,
300                                     bShowNumbers ? start : -1, bShowNumbers ? end - 1 : -1);
301                 }
302                 for (j = start; j < end; j++)
303                 {
304                     if (j > start)
305                     {
306                         size += fprintf(fp, ", ");
307                     }
308                     if ((size) > (USE_WIDTH))
309                     {
310                         fprintf(fp, "\n");
311                         size = pr_indent(fp, indent + INDENT);
312                     }
313                     size += fprintf(fp, "%d", block->a[j]);
314                 }
315                 fprintf(fp, "}\n");
316                 start = end;
317             }
318         }
319         if ((end != block->nra) || (!ok))
320         {
321             pr_indent(fp, indent);
322             fprintf(fp, "tables inconsistent, dumping complete tables:\n");
323             low_pr_blocka(fp, indent, title, block, bShowNumbers);
324         }
325     }
326 }
327
328 void copy_block(const t_block* src, t_block* dst)
329 {
330     dst->nr = src->nr;
331     /* Workaround for inconsistent handling of nalloc_index in
332      * other parts of the code. Often nalloc_index and nalloc_a
333      * are not set.
334      */
335     dst->nalloc_index = std::max(src->nalloc_index, dst->nr + 1);
336     snew(dst->index, dst->nalloc_index);
337     for (int i = 0; i < dst->nr + 1; ++i)
338     {
339         dst->index[i] = src->index[i];
340     }
341 }