8784374c3d96bc7a2e78b272f70f6e682eb36277
[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, 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_blocka(t_blocka *block)
76 {
77     block->nr           = 0;
78     block->nra          = 0;
79     block->nalloc_index = 1;
80     snew(block->index, block->nalloc_index);
81     block->index[0]     = 0;
82     block->nalloc_a     = 0;
83     block->a            = nullptr;
84 }
85
86 t_blocka *new_blocka(void)
87 {
88     t_blocka *block;
89
90     snew(block, 1);
91     snew(block->index, 1);
92
93     return block;
94 }
95
96 void done_block(t_block *block)
97 {
98     block->nr           = 0;
99     sfree(block->index);
100     block->index        = nullptr;
101     block->nalloc_index = 0;
102 }
103
104 void done_blocka(t_blocka *block)
105 {
106     block->nr    = 0;
107     block->nra   = 0;
108     sfree(block->index);
109     sfree(block->a);
110     block->index        = nullptr;
111     block->a            = nullptr;
112     block->nalloc_index = 0;
113     block->nalloc_a     = 0;
114 }
115
116 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup)
117 {
118     if (bOneIndexGroup)
119     {
120         grp->nalloc_index = 2;
121         snew(grp->index, grp->nalloc_index);
122         grp->index[0] = 0;
123         grp->index[1] = natom;
124         grp->nr       = 1;
125     }
126     else
127     {
128         grp->nalloc_index = natom+1;
129         snew(grp->index, grp->nalloc_index);
130         snew(grp->index, natom+1);
131         for (int i = 0; i <= natom; ++i)
132         {
133             grp->index[i] = i;
134         }
135         grp->nr = natom;
136     }
137 }
138
139 void stupid_fill_blocka(t_blocka *grp, int natom)
140 {
141     grp->nalloc_a = natom;
142     snew(grp->a, grp->nalloc_a);
143     for (int i = 0; i < natom; ++i)
144     {
145         grp->a[i] = i;
146     }
147     grp->nra = natom;
148
149     grp->nalloc_index = natom + 1;
150     snew(grp->index, grp->nalloc_index);
151     for (int i = 0; i <= natom; ++i)
152     {
153         grp->index[i] = i;
154     }
155     grp->nr = natom;
156 }
157
158 void copy_blocka(const t_blocka *src, t_blocka *dest)
159 {
160     dest->nr           = src->nr;
161     /* Workaround for inconsistent handling of nalloc_index in
162      * other parts of the code. Often nalloc_index and nalloc_a
163      * are not set.
164      */
165     dest->nalloc_index = std::max(src->nalloc_index, dest->nr + 1);
166     snew(dest->index, dest->nalloc_index);
167     for (int i = 0; i < dest->nr+1; ++i)
168     {
169         dest->index[i] = src->index[i];
170     }
171     dest->nra      = src->nra;
172     /* See above. */
173     dest->nalloc_a = std::max(src->nalloc_a, dest->nra);
174     snew(dest->a, dest->nalloc_a);
175     for (int i = 0; i < dest->nra; ++i)
176     {
177         dest->a[i] = src->a[i];
178     }
179 }
180
181 static int pr_block_title(FILE *fp, int indent, const char *title, const t_block *block)
182 {
183     if (available(fp, block, indent, title))
184     {
185         indent = pr_title(fp, indent, title);
186         pr_indent(fp, indent);
187         fprintf(fp, "nr=%d\n", block->nr);
188     }
189     return indent;
190 }
191
192 static int pr_blocka_title(FILE *fp, int indent, const char *title, const t_blocka *block)
193 {
194     if (available(fp, block, indent, title))
195     {
196         indent = pr_title(fp, indent, title);
197         pr_indent(fp, indent);
198         fprintf(fp, "nr=%d\n", block->nr);
199         pr_indent(fp, indent);
200         fprintf(fp, "nra=%d\n", block->nra);
201     }
202     return indent;
203 }
204
205 static void low_pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers)
206 {
207     int i;
208
209     if (available(fp, block, indent, title))
210     {
211         indent = pr_blocka_title(fp, indent, title, block);
212         for (i = 0; i <= block->nr; i++)
213         {
214             pr_indent(fp, indent+INDENT);
215             fprintf(fp, "%s->index[%d]=%d\n",
216                     title, bShowNumbers ? i : -1, block->index[i]);
217         }
218         for (i = 0; i < block->nra; i++)
219         {
220             pr_indent(fp, indent+INDENT);
221             fprintf(fp, "%s->a[%d]=%d\n",
222                     title, bShowNumbers ? i : -1, block->a[i]);
223         }
224     }
225 }
226
227 void pr_block(FILE *fp, int indent, const char *title, const t_block *block, gmx_bool bShowNumbers)
228 {
229     int i, start;
230
231     if (available(fp, block, indent, title))
232     {
233         indent = pr_block_title(fp, indent, title, block);
234         start  = 0;
235         if (block->index[start] != 0)
236         {
237             fprintf(fp, "block->index[%d] should be 0\n", start);
238         }
239         else
240         {
241             for (i = 0; i < block->nr; i++)
242             {
243                 int end  = block->index[i+1];
244                 pr_indent(fp, indent);
245                 if (end <= start)
246                 {
247                     fprintf(fp, "%s[%d]={}\n", title, i);
248                 }
249                 else
250                 {
251                     fprintf(fp, "%s[%d]={%d..%d}\n",
252                             title, bShowNumbers ? i : -1,
253                             bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
254                 }
255                 start = end;
256             }
257         }
258     }
259 }
260
261 void pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers)
262 {
263     int i, j, ok, size, start, end;
264
265     if (available(fp, block, indent, title))
266     {
267         indent = pr_blocka_title(fp, indent, title, block);
268         start  = 0;
269         end    = start;
270         if ((ok = (block->index[start] == 0)) == 0)
271         {
272             fprintf(fp, "block->index[%d] should be 0\n", start);
273         }
274         else
275         {
276             for (i = 0; i < block->nr; i++)
277             {
278                 end  = block->index[i+1];
279                 size = pr_indent(fp, indent);
280                 if (end <= start)
281                 {
282                     size += fprintf(fp, "%s[%d]={", title, i);
283                 }
284                 else
285                 {
286                     size += fprintf(fp, "%s[%d][%d..%d]={",
287                                     title, bShowNumbers ? i : -1,
288                                     bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
289                 }
290                 for (j = start; j < end; j++)
291                 {
292                     if (j > start)
293                     {
294                         size += fprintf(fp, ", ");
295                     }
296                     if ((size) > (USE_WIDTH))
297                     {
298                         fprintf(fp, "\n");
299                         size = pr_indent(fp, indent+INDENT);
300                     }
301                     size += fprintf(fp, "%d", block->a[j]);
302                 }
303                 fprintf(fp, "}\n");
304                 start = end;
305             }
306         }
307         if ((end != block->nra) || (!ok))
308         {
309             pr_indent(fp, indent);
310             fprintf(fp, "tables inconsistent, dumping complete tables:\n");
311             low_pr_blocka(fp, indent, title, block, bShowNumbers);
312         }
313     }
314 }