2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #include "sparsematrix.h"
43 #include "gromacs/legacyheaders/smalloc.h"
49 gmx_sparsematrix_init(int nrow)
52 gmx_sparsematrix_t *A;
58 snew(A->nalloc, nrow);
61 for (i = 0; i < nrow; i++)
73 gmx_sparsematrix_destroy(gmx_sparsematrix_t * A)
77 /* Release each row */
78 for (i = 0; i < A->nrow; i++)
80 if (A->data[i] != NULL)
85 /* Release the rowdata arrays */
89 /* Release matrix structure itself */
96 gmx_sparsematrix_print(FILE * stream,
97 gmx_sparsematrix_t * A)
101 for (i = 0; i < A->nrow; i++)
103 if (A->ndata[i] == 0)
105 for (j = 0; j < A->nrow; j++)
107 fprintf(stream, " %6.3f", 0.0);
114 for (j = 0; j < A->ndata[i]; j++)
116 while (k++ < A->data[i][j].col)
118 fprintf(stream, " %6.3f", 0.0);
120 fprintf(stream, " %6.3f", A->data[i][j].value);
122 while (k++ < A->nrow)
124 fprintf(stream, " %6.3f", 0.0);
127 fprintf(stream, "\n");
134 gmx_sparsematrix_value(gmx_sparsematrix_t * A,
138 gmx_bool found = FALSE;
142 assert(row < A->nrow);
146 /* Find previous value */
147 for (i = 0; i < A->ndata[row] && (found == FALSE); i++)
149 if (A->data[row][i].col == col)
152 value = A->data[row][i].value;
156 /* value=0 if we didn't find any match */
163 gmx_sparsematrix_increment_value(gmx_sparsematrix_t * A,
168 gmx_bool found = FALSE;
171 assert(row < A->nrow);
173 /* Try to find a previous entry with this row/col */
174 for (i = 0; i < A->ndata[row] && !found; i++)
176 if (A->data[row][i].col == col)
179 A->data[row][i].value += difference;
183 /* Add a new entry if nothing was found */
186 /* add the value at the end of the row */
187 if (A->ndata[row] == A->nalloc[row])
189 A->nalloc[row] += 100;
190 if (A->data[row] == NULL)
192 snew(A->data[row], A->nalloc[row]);
196 srenew(A->data[row], A->nalloc[row]);
199 A->data[row][A->ndata[row]].col = col;
200 /* Previous value was 0.0 */
201 A->data[row][A->ndata[row]].value = difference;
207 /* Routine to compare column values of two entries, used for quicksort of each row.
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.
213 compare_columns(const void *v1, const void *v2)
215 int c1 = ((gmx_sparsematrix_entry_t *)v1)->col;
216 int c2 = ((gmx_sparsematrix_entry_t *)v2)->col;
234 gmx_sparsematrix_compress(gmx_sparsematrix_t * A)
238 for (i = 0; i < A->nrow; i++)
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)
246 /* Go through values on this row and look for more zero elements */
247 for (j = 0; j < A->ndata[i]; j++)
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)
252 if (A->data[i][j].value == 0)
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;
259 /* Only non-zero elements remaining on this row. Sort them after column index */
260 qsort((void *)(A->data[i]),
262 sizeof(gmx_sparsematrix_entry_t),
269 gmx_sparsematrix_vector_multiply(gmx_sparsematrix_t * A,
275 gmx_sparsematrix_entry_t * data; /* pointer to simplify data access */
277 for (i = 0; i < A->nrow; i++)
282 if (A->compressed_symmetric)
284 for (i = 0; i < A->nrow; i++)
290 for (k = 0; k < A->ndata[i]; k++)
305 /* not compressed symmetric storage */
306 for (i = 0; i < A->nrow; i++)
312 for (k = 0; k < A->ndata[i]; k++)