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