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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "sparsematrix.h"
54 gmx_sparsematrix_init(int nrow)
57 gmx_sparsematrix_t *A;
78 gmx_sparsematrix_destroy(gmx_sparsematrix_t * A)
82 /* Release each row */
83 for(i=0;i<A->nrow;i++)
88 /* Release the rowdata arrays */
92 /* Release matrix structure itself */
99 gmx_sparsematrix_print(FILE * stream,
100 gmx_sparsematrix_t * A)
104 for(i=0;i<A->nrow;i++)
108 for(j=0;j<A->nrow;j++)
109 fprintf(stream," %6.3f",0.0);
115 for(j=0;j<A->ndata[i];j++)
117 while(k++<A->data[i][j].col)
118 fprintf(stream," %6.3f",0.0);
119 fprintf(stream," %6.3f",A->data[i][j].value);
122 fprintf(stream," %6.3f",0.0);
124 fprintf(stream,"\n");
131 gmx_sparsematrix_value(gmx_sparsematrix_t * A,
135 gmx_bool found = FALSE;
143 /* Find previous value */
144 for(i=0;i<A->ndata[row] && (found==FALSE);i++)
146 if(A->data[row][i].col==col)
149 value = A->data[row][i].value;
153 /* value=0 if we didn't find any match */
160 gmx_sparsematrix_increment_value(gmx_sparsematrix_t * A,
165 gmx_bool found = FALSE;
170 /* Try to find a previous entry with this row/col */
171 for(i=0;i<A->ndata[row] && !found;i++)
173 if(A->data[row][i].col==col)
176 A->data[row][i].value += difference;
180 /* Add a new entry if nothing was found */
183 /* add the value at the end of the row */
184 if(A->ndata[row] == A->nalloc[row])
186 A->nalloc[row] += 100;
187 if(A->data[row]==NULL)
189 snew(A->data[row],A->nalloc[row]);
193 srenew(A->data[row],A->nalloc[row]);
196 A->data[row][A->ndata[row]].col = col;
197 /* Previous value was 0.0 */
198 A->data[row][A->ndata[row]].value = difference;
204 /* Routine to compare column values of two entries, used for quicksort of each row.
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.
210 compare_columns(const void *v1, const void *v2)
212 int c1 = ((gmx_sparsematrix_entry_t *)v1)->col;
213 int c2 = ((gmx_sparsematrix_entry_t *)v2)->col;
225 gmx_sparsematrix_compress(gmx_sparsematrix_t * A)
229 for (i=0;i<A->nrow;i++)
231 /* Remove last value on this row while it is zero */
232 while(A->ndata[i]>0 && A->data[i][A->ndata[i]-1].value==0)
235 /* Go through values on this row and look for more zero elements */
236 for(j=0;j<A->ndata[i];j++)
238 /* If this element was zero, exchange it with the last non-zero
239 * element on the row (yes, this will invalidate the sort order)
241 if(A->data[i][j].value==0)
243 A->data[i][j].value = A->data[i][A->ndata[i]-1].value;
244 A->data[i][j].col = A->data[i][A->ndata[i]-1].col;
248 /* Only non-zero elements remaining on this row. Sort them after column index */
249 qsort((void *)(A->data[i]),
251 sizeof(gmx_sparsematrix_entry_t),
258 gmx_sparsematrix_vector_multiply(gmx_sparsematrix_t * A,
264 gmx_sparsematrix_entry_t * data; /* pointer to simplify data access */
266 for (i = 0; i < A->nrow; i ++)
269 if(A->compressed_symmetric)
271 for (i = 0; i < A->nrow; i ++)
277 for (k=0;k<A->ndata[i];k++)
290 /* not compressed symmetric storage */
291 for (i = 0; i < A->nrow; i ++)
297 for (k=0;k<A->ndata[i];k++)