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