Merge branch release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / linearalgebra / sparsematrix.c
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) 2012,2014, 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 "sparsematrix.h"
38
39 #include <assert.h>
40 #include <stdlib.h>
41 #include <stdio.h>
42
43 #include "gromacs/utility/smalloc.h"
44
45 gmx_sparsematrix_t *
46 gmx_sparsematrix_init(int                    nrow)
47 {
48     int                 i;
49     gmx_sparsematrix_t *A;
50
51     snew(A, 1);
52
53     A->nrow = nrow;
54     snew(A->ndata, nrow);
55     snew(A->nalloc, nrow);
56     snew(A->data, nrow);
57
58     for (i = 0; i < nrow; i++)
59     {
60         A->ndata[i]    = 0;
61         A->nalloc[i]   = 0;
62         A->data[i]     = NULL;
63     }
64     return A;
65 }
66
67
68
69 void
70 gmx_sparsematrix_destroy(gmx_sparsematrix_t *   A)
71 {
72     int i;
73
74     /* Release each row */
75     for (i = 0; i < A->nrow; i++)
76     {
77         if (A->data[i] != NULL)
78         {
79             sfree(A->data[i]);
80         }
81     }
82     /* Release the rowdata arrays */
83     sfree(A->ndata);
84     sfree(A->nalloc);
85     sfree(A->data);
86     /* Release matrix structure itself */
87     sfree(A);
88 }
89
90
91
92 void
93 gmx_sparsematrix_print(FILE *                 stream,
94                        gmx_sparsematrix_t *   A)
95 {
96     int i, j, k;
97
98     for (i = 0; i < A->nrow; i++)
99     {
100         if (A->ndata[i] == 0)
101         {
102             for (j = 0; j < A->nrow; j++)
103             {
104                 fprintf(stream, " %6.3f", 0.0);
105             }
106         }
107         else
108         {
109             k = 0;
110             j = 0;
111             for (j = 0; j < A->ndata[i]; j++)
112             {
113                 while (k++ < A->data[i][j].col)
114                 {
115                     fprintf(stream, " %6.3f", 0.0);
116                 }
117                 fprintf(stream, " %6.3f", A->data[i][j].value);
118             }
119             while (k++ < A->nrow)
120             {
121                 fprintf(stream, " %6.3f", 0.0);
122             }
123         }
124         fprintf(stream, "\n");
125     }
126
127 }
128
129
130 real
131 gmx_sparsematrix_value(gmx_sparsematrix_t *    A,
132                        int                     row,
133                        int                     col)
134 {
135     gmx_bool found  = FALSE;
136     int      i;
137     real     value;
138
139     assert(row < A->nrow);
140
141     value = 0;
142
143     /* Find previous value */
144     for (i = 0; i < A->ndata[row] && (found == FALSE); i++)
145     {
146         if (A->data[row][i].col == col)
147         {
148             found = TRUE;
149             value = A->data[row][i].value;
150         }
151     }
152
153     /* value=0 if we didn't find any match */
154     return value;
155 }
156
157
158
159 void
160 gmx_sparsematrix_increment_value(gmx_sparsematrix_t *    A,
161                                  int                     row,
162                                  int                     col,
163                                  real                    difference)
164 {
165     gmx_bool found  = FALSE;
166     int      i;
167
168     assert(row < A->nrow);
169
170     /* Try to find a previous entry with this row/col */
171     for (i = 0; i < A->ndata[row] && !found; i++)
172     {
173         if (A->data[row][i].col == col)
174         {
175             found                  = TRUE;
176             A->data[row][i].value += difference;
177         }
178     }
179
180     /* Add a new entry if nothing was found */
181     if (!found)
182     {
183         /* add the value at the end of the row */
184         if (A->ndata[row] == A->nalloc[row])
185         {
186             A->nalloc[row] += 100;
187             if (A->data[row] == NULL)
188             {
189                 snew(A->data[row], A->nalloc[row]);
190             }
191             else
192             {
193                 srenew(A->data[row], A->nalloc[row]);
194             }
195         }
196         A->data[row][A->ndata[row]].col   = col;
197         /* Previous value was 0.0 */
198         A->data[row][A->ndata[row]].value = difference;
199         A->ndata[row]++;
200     }
201 }
202
203
204 /* Routine to compare column values of two entries, used for quicksort of each row.
205  *
206  * The data entries to compare are of the type gmx_sparsematrix_entry_t, but quicksort
207  * uses void pointers as arguments, so we cast them back internally.
208  */
209 static int
210 compare_columns(const void *v1, const void *v2)
211 {
212     int c1 = ((gmx_sparsematrix_entry_t *)v1)->col;
213     int c2 = ((gmx_sparsematrix_entry_t *)v2)->col;
214
215     if (c1 < c2)
216     {
217         return -1;
218     }
219     else if (c1 > c2)
220     {
221         return 1;
222     }
223     else
224     {
225         return 0;
226     }
227 }
228
229
230 void
231 gmx_sparsematrix_compress(gmx_sparsematrix_t *    A)
232 {
233     int i, j;
234
235     for (i = 0; i < A->nrow; i++)
236     {
237         /* Remove last value on this row while it is zero */
238         while (A->ndata[i] > 0 && A->data[i][A->ndata[i]-1].value == 0)
239         {
240             A->ndata[i]--;
241         }
242
243         /* Go through values on this row and look for more zero elements */
244         for (j = 0; j < A->ndata[i]; j++)
245         {
246             /* If this element was zero, exchange it with the last non-zero
247              * element on the row (yes, this will invalidate the sort order)
248              */
249             if (A->data[i][j].value == 0)
250             {
251                 A->data[i][j].value = A->data[i][A->ndata[i]-1].value;
252                 A->data[i][j].col   = A->data[i][A->ndata[i]-1].col;
253                 A->ndata[i]--;
254             }
255         }
256         /* Only non-zero elements remaining on this row. Sort them after column index */
257         qsort((void *)(A->data[i]),
258               A->ndata[i],
259               sizeof(gmx_sparsematrix_entry_t),
260               compare_columns);
261     }
262 }
263
264
265 void
266 gmx_sparsematrix_vector_multiply(gmx_sparsematrix_t *    A,
267                                  real *                  x,
268                                  real *                  y)
269 {
270     real                        s, v, xi;
271     int                         i, j, k;
272     gmx_sparsematrix_entry_t *  data; /* pointer to simplify data access */
273
274     for (i = 0; i < A->nrow; i++)
275     {
276         y[i] = 0;
277     }
278
279     if (A->compressed_symmetric)
280     {
281         for (i = 0; i < A->nrow; i++)
282         {
283             xi   = x[i];
284             s    = 0.0;
285             data = A->data[i];
286
287             for (k = 0; k < A->ndata[i]; k++)
288             {
289                 j  = data[k].col;
290                 v  = data[k].value;
291                 s += v * x[j];
292                 if (i != j)
293                 {
294                     y[j] += v * xi;
295                 }
296             }
297             y[i] += s;
298         }
299     }
300     else
301     {
302         /* not compressed symmetric storage */
303         for (i = 0; i < A->nrow; i++)
304         {
305             xi   = x[i];
306             s    = 0.0;
307             data = A->data[i];
308
309             for (k = 0; k < A->ndata[i]; k++)
310             {
311                 j  = data[k].col;
312                 v  = data[k].value;
313                 s += v * x[j];
314             }
315             y[i] += s;
316         }
317     }
318 }