3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
45 #include "sparsematrix.h"
51 gmx_sparsematrix_init(int nrow
)
54 gmx_sparsematrix_t
*A
;
75 gmx_sparsematrix_destroy(gmx_sparsematrix_t
* A
)
79 /* Release each row */
80 for(i
=0;i
<A
->nrow
;i
++)
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
++)
105 for(j
=0;j
<A
->nrow
;j
++)
106 fprintf(stream
," %6.3f",0.0);
112 for(j
=0;j
<A
->ndata
[i
];j
++)
114 while(k
++<A
->data
[i
][j
].col
)
115 fprintf(stream
," %6.3f",0.0);
116 fprintf(stream
," %6.3f",A
->data
[i
][j
].value
);
119 fprintf(stream
," %6.3f",0.0);
121 fprintf(stream
,"\n");
128 gmx_sparsematrix_value(gmx_sparsematrix_t
* A
,
132 gmx_bool found
= FALSE
;
140 /* Find previous value */
141 for(i
=0;i
<A
->ndata
[row
] && (found
==FALSE
);i
++)
143 if(A
->data
[row
][i
].col
==col
)
146 value
= A
->data
[row
][i
].value
;
150 /* value=0 if we didn't find any match */
157 gmx_sparsematrix_increment_value(gmx_sparsematrix_t
* A
,
162 gmx_bool found
= FALSE
;
167 /* Try to find a previous entry with this row/col */
168 for(i
=0;i
<A
->ndata
[row
] && !found
;i
++)
170 if(A
->data
[row
][i
].col
==col
)
173 A
->data
[row
][i
].value
+= difference
;
177 /* Add a new entry if nothing was found */
180 /* add the value at the end of the row */
181 if(A
->ndata
[row
] == A
->nalloc
[row
])
183 A
->nalloc
[row
] += 100;
184 if(A
->data
[row
]==NULL
)
186 snew(A
->data
[row
],A
->nalloc
[row
]);
190 srenew(A
->data
[row
],A
->nalloc
[row
]);
193 A
->data
[row
][A
->ndata
[row
]].col
= col
;
194 /* Previous value was 0.0 */
195 A
->data
[row
][A
->ndata
[row
]].value
= difference
;
201 /* Routine to compare column values of two entries, used for quicksort of each row.
203 * The data entries to compare are of the type gmx_sparsematrix_entry_t, but quicksort
204 * uses void pointers as arguments, so we cast them back internally.
207 compare_columns(const void *v1
, const void *v2
)
209 int c1
= ((gmx_sparsematrix_entry_t
*)v1
)->col
;
210 int c2
= ((gmx_sparsematrix_entry_t
*)v2
)->col
;
222 gmx_sparsematrix_compress(gmx_sparsematrix_t
* A
)
226 for (i
=0;i
<A
->nrow
;i
++)
228 /* Remove last value on this row while it is zero */
229 while(A
->ndata
[i
]>0 && A
->data
[i
][A
->ndata
[i
]-1].value
==0)
232 /* Go through values on this row and look for more zero elements */
233 for(j
=0;j
<A
->ndata
[i
];j
++)
235 /* If this element was zero, exchange it with the last non-zero
236 * element on the row (yes, this will invalidate the sort order)
238 if(A
->data
[i
][j
].value
==0)
240 A
->data
[i
][j
].value
= A
->data
[i
][A
->ndata
[i
]-1].value
;
241 A
->data
[i
][j
].col
= A
->data
[i
][A
->ndata
[i
]-1].col
;
245 /* Only non-zero elements remaining on this row. Sort them after column index */
246 qsort((void *)(A
->data
[i
]),
248 sizeof(gmx_sparsematrix_entry_t
),
255 gmx_sparsematrix_vector_multiply(gmx_sparsematrix_t
* A
,
261 gmx_sparsematrix_entry_t
* data
; /* pointer to simplify data access */
263 for (i
= 0; i
< A
->nrow
; i
++)
266 if(A
->compressed_symmetric
)
268 for (i
= 0; i
< A
->nrow
; i
++)
274 for (k
=0;k
<A
->ndata
[i
];k
++)
287 /* not compressed symmetric storage */
288 for (i
= 0; i
< A
->nrow
; i
++)
294 for (k
=0;k
<A
->ndata
[i
];k
++)