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