alpar@1: /* glpspm.c */ alpar@1: alpar@1: /*********************************************************************** alpar@1: * This code is part of GLPK (GNU Linear Programming Kit). alpar@1: * alpar@1: * Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, alpar@1: * 2009, 2010 Andrew Makhorin, Department for Applied Informatics, alpar@1: * Moscow Aviation Institute, Moscow, Russia. All rights reserved. alpar@1: * E-mail: . alpar@1: * alpar@1: * GLPK is free software: you can redistribute it and/or modify it alpar@1: * under the terms of the GNU General Public License as published by alpar@1: * the Free Software Foundation, either version 3 of the License, or alpar@1: * (at your option) any later version. alpar@1: * alpar@1: * GLPK is distributed in the hope that it will be useful, but WITHOUT alpar@1: * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY alpar@1: * or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public alpar@1: * License for more details. alpar@1: * alpar@1: * You should have received a copy of the GNU General Public License alpar@1: * along with GLPK. If not, see . alpar@1: ***********************************************************************/ alpar@1: alpar@1: #include "glphbm.h" alpar@1: #include "glprgr.h" alpar@1: #include "glpspm.h" alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_create_mat - create general sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_create_mat(int m, int n); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_create_mat creates a general sparse matrix having alpar@1: * m rows and n columns. Being created the matrix is zero (empty), i.e. alpar@1: * has no elements. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine returns a pointer to the matrix created. */ alpar@1: alpar@1: SPM *spm_create_mat(int m, int n) alpar@1: { SPM *A; alpar@1: xassert(0 <= m && m < INT_MAX); alpar@1: xassert(0 <= n && n < INT_MAX); alpar@1: A = xmalloc(sizeof(SPM)); alpar@1: A->m = m; alpar@1: A->n = n; alpar@1: if (m == 0 || n == 0) alpar@1: { A->pool = NULL; alpar@1: A->row = NULL; alpar@1: A->col = NULL; alpar@1: } alpar@1: else alpar@1: { int i, j; alpar@1: A->pool = dmp_create_pool(); alpar@1: A->row = xcalloc(1+m, sizeof(SPME *)); alpar@1: for (i = 1; i <= m; i++) A->row[i] = NULL; alpar@1: A->col = xcalloc(1+n, sizeof(SPME *)); alpar@1: for (j = 1; j <= n; j++) A->col[j] = NULL; alpar@1: } alpar@1: return A; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_new_elem - add new element to sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPME *spm_new_elem(SPM *A, int i, int j, double val); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_new_elem adds a new element to the specified sparse alpar@1: * matrix. Parameters i, j, and val specify the row number, the column alpar@1: * number, and a numerical value of the element, respectively. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine returns a pointer to the new element added. */ alpar@1: alpar@1: SPME *spm_new_elem(SPM *A, int i, int j, double val) alpar@1: { SPME *e; alpar@1: xassert(1 <= i && i <= A->m); alpar@1: xassert(1 <= j && j <= A->n); alpar@1: e = dmp_get_atom(A->pool, sizeof(SPME)); alpar@1: e->i = i; alpar@1: e->j = j; alpar@1: e->val = val; alpar@1: e->r_prev = NULL; alpar@1: e->r_next = A->row[i]; alpar@1: if (e->r_next != NULL) e->r_next->r_prev = e; alpar@1: e->c_prev = NULL; alpar@1: e->c_next = A->col[j]; alpar@1: if (e->c_next != NULL) e->c_next->c_prev = e; alpar@1: A->row[i] = A->col[j] = e; alpar@1: return e; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_delete_mat - delete general sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * void spm_delete_mat(SPM *A); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine deletes the specified general sparse matrix freeing all alpar@1: * the memory allocated to this object. */ alpar@1: alpar@1: void spm_delete_mat(SPM *A) alpar@1: { /* delete sparse matrix */ alpar@1: if (A->pool != NULL) dmp_delete_pool(A->pool); alpar@1: if (A->row != NULL) xfree(A->row); alpar@1: if (A->col != NULL) xfree(A->col); alpar@1: xfree(A); alpar@1: return; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_test_mat_e - create test sparse matrix of E(n,c) class alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_test_mat_e(int n, int c); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_test_mat_e creates a test sparse matrix of E(n,c) alpar@1: * class as described in the book: Ole 0sterby, Zahari Zlatev. Direct alpar@1: * Methods for Sparse Matrices. Springer-Verlag, 1983. alpar@1: * alpar@1: * Matrix of E(n,c) class is a symmetric positive definite matrix of alpar@1: * the order n. It has the number 4 on its main diagonal and the number alpar@1: * -1 on its four co-diagonals, two of which are neighbour to the main alpar@1: * diagonal and two others are shifted from the main diagonal on the alpar@1: * distance c. alpar@1: * alpar@1: * It is necessary that n >= 3 and 2 <= c <= n-1. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine returns a pointer to the matrix created. */ alpar@1: alpar@1: SPM *spm_test_mat_e(int n, int c) alpar@1: { SPM *A; alpar@1: int i; alpar@1: xassert(n >= 3 && 2 <= c && c <= n-1); alpar@1: A = spm_create_mat(n, n); alpar@1: for (i = 1; i <= n; i++) alpar@1: spm_new_elem(A, i, i, 4.0); alpar@1: for (i = 1; i <= n-1; i++) alpar@1: { spm_new_elem(A, i, i+1, -1.0); alpar@1: spm_new_elem(A, i+1, i, -1.0); alpar@1: } alpar@1: for (i = 1; i <= n-c; i++) alpar@1: { spm_new_elem(A, i, i+c, -1.0); alpar@1: spm_new_elem(A, i+c, i, -1.0); alpar@1: } alpar@1: return A; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_test_mat_d - create test sparse matrix of D(n,c) class alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_test_mat_d(int n, int c); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_test_mat_d creates a test sparse matrix of D(n,c) alpar@1: * class as described in the book: Ole 0sterby, Zahari Zlatev. Direct alpar@1: * Methods for Sparse Matrices. Springer-Verlag, 1983. alpar@1: * alpar@1: * Matrix of D(n,c) class is a non-singular matrix of the order n. It alpar@1: * has unity main diagonal, three co-diagonals above the main diagonal alpar@1: * on the distance c, which are cyclically continued below the main alpar@1: * diagonal, and a triangle block of the size 10x10 in the upper right alpar@1: * corner. alpar@1: * alpar@1: * It is necessary that n >= 14 and 1 <= c <= n-13. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine returns a pointer to the matrix created. */ alpar@1: alpar@1: SPM *spm_test_mat_d(int n, int c) alpar@1: { SPM *A; alpar@1: int i, j; alpar@1: xassert(n >= 14 && 1 <= c && c <= n-13); alpar@1: A = spm_create_mat(n, n); alpar@1: for (i = 1; i <= n; i++) alpar@1: spm_new_elem(A, i, i, 1.0); alpar@1: for (i = 1; i <= n-c; i++) alpar@1: spm_new_elem(A, i, i+c, (double)(i+1)); alpar@1: for (i = n-c+1; i <= n; i++) alpar@1: spm_new_elem(A, i, i-n+c, (double)(i+1)); alpar@1: for (i = 1; i <= n-c-1; i++) alpar@1: spm_new_elem(A, i, i+c+1, (double)(-i)); alpar@1: for (i = n-c; i <= n; i++) alpar@1: spm_new_elem(A, i, i-n+c+1, (double)(-i)); alpar@1: for (i = 1; i <= n-c-2; i++) alpar@1: spm_new_elem(A, i, i+c+2, 16.0); alpar@1: for (i = n-c-1; i <= n; i++) alpar@1: spm_new_elem(A, i, i-n+c+2, 16.0); alpar@1: for (j = 1; j <= 10; j++) alpar@1: for (i = 1; i <= 11-j; i++) alpar@1: spm_new_elem(A, i, n-11+i+j, 100.0 * (double)j); alpar@1: return A; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_show_mat - write sparse matrix pattern in BMP file format alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * int spm_show_mat(const SPM *A, const char *fname); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_show_mat writes pattern of the specified sparse alpar@1: * matrix in uncompressed BMP file format (Windows bitmap) to a binary alpar@1: * file whose name is specified by the character string fname. alpar@1: * alpar@1: * Each pixel corresponds to one matrix element. The pixel colors have alpar@1: * the following meaning: alpar@1: * alpar@1: * Black structurally zero element alpar@1: * White positive element alpar@1: * Cyan negative element alpar@1: * Green zero element alpar@1: * Red duplicate element alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * If no error occured, the routine returns zero. Otherwise, it prints alpar@1: * an appropriate error message and returns non-zero. */ alpar@1: alpar@1: int spm_show_mat(const SPM *A, const char *fname) alpar@1: { int m = A->m; alpar@1: int n = A->n; alpar@1: int i, j, k, ret; alpar@1: char *map; alpar@1: xprintf("spm_show_mat: writing matrix pattern to `%s'...\n", alpar@1: fname); alpar@1: xassert(1 <= m && m <= 32767); alpar@1: xassert(1 <= n && n <= 32767); alpar@1: map = xmalloc(m * n); alpar@1: memset(map, 0x08, m * n); alpar@1: for (i = 1; i <= m; i++) alpar@1: { SPME *e; alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: { j = e->j; alpar@1: xassert(1 <= j && j <= n); alpar@1: k = n * (i - 1) + (j - 1); alpar@1: if (map[k] != 0x08) alpar@1: map[k] = 0x0C; alpar@1: else if (e->val > 0.0) alpar@1: map[k] = 0x0F; alpar@1: else if (e->val < 0.0) alpar@1: map[k] = 0x0B; alpar@1: else alpar@1: map[k] = 0x0A; alpar@1: } alpar@1: } alpar@1: ret = rgr_write_bmp16(fname, m, n, map); alpar@1: xfree(map); alpar@1: return ret; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_read_hbm - read sparse matrix in Harwell-Boeing format alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_read_hbm(const char *fname); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_read_hbm reads a sparse matrix in the Harwell-Boeing alpar@1: * format from a text file whose name is the character string fname. alpar@1: * alpar@1: * Detailed description of the Harwell-Boeing format recognised by this alpar@1: * routine can be found in the following report: alpar@1: * alpar@1: * I.S.Duff, R.G.Grimes, J.G.Lewis. User's Guide for the Harwell-Boeing alpar@1: * Sparse Matrix Collection (Release I), TR/PA/92/86, October 1992. alpar@1: * alpar@1: * NOTE alpar@1: * alpar@1: * The routine spm_read_hbm reads the matrix "as is", due to which zero alpar@1: * and/or duplicate elements can appear in the matrix. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * If no error occured, the routine returns a pointer to the matrix alpar@1: * created. Otherwise, the routine prints an appropriate error message alpar@1: * and returns NULL. */ alpar@1: alpar@1: SPM *spm_read_hbm(const char *fname) alpar@1: { SPM *A = NULL; alpar@1: HBM *hbm; alpar@1: int nrow, ncol, nnzero, i, j, beg, end, ptr, *colptr, *rowind; alpar@1: double val, *values; alpar@1: char *mxtype; alpar@1: hbm = hbm_read_mat(fname); alpar@1: if (hbm == NULL) alpar@1: { xprintf("spm_read_hbm: unable to read matrix\n"); alpar@1: goto fini; alpar@1: } alpar@1: mxtype = hbm->mxtype; alpar@1: nrow = hbm->nrow; alpar@1: ncol = hbm->ncol; alpar@1: nnzero = hbm->nnzero; alpar@1: colptr = hbm->colptr; alpar@1: rowind = hbm->rowind; alpar@1: values = hbm->values; alpar@1: if (!(strcmp(mxtype, "RSA") == 0 || strcmp(mxtype, "PSA") == 0 || alpar@1: strcmp(mxtype, "RUA") == 0 || strcmp(mxtype, "PUA") == 0 || alpar@1: strcmp(mxtype, "RRA") == 0 || strcmp(mxtype, "PRA") == 0)) alpar@1: { xprintf("spm_read_hbm: matrix type `%s' not supported\n", alpar@1: mxtype); alpar@1: goto fini; alpar@1: } alpar@1: A = spm_create_mat(nrow, ncol); alpar@1: if (mxtype[1] == 'S' || mxtype[1] == 'U') alpar@1: xassert(nrow == ncol); alpar@1: for (j = 1; j <= ncol; j++) alpar@1: { beg = colptr[j]; alpar@1: end = colptr[j+1]; alpar@1: xassert(1 <= beg && beg <= end && end <= nnzero + 1); alpar@1: for (ptr = beg; ptr < end; ptr++) alpar@1: { i = rowind[ptr]; alpar@1: xassert(1 <= i && i <= nrow); alpar@1: if (mxtype[0] == 'R') alpar@1: val = values[ptr]; alpar@1: else alpar@1: val = 1.0; alpar@1: spm_new_elem(A, i, j, val); alpar@1: if (mxtype[1] == 'S' && i != j) alpar@1: spm_new_elem(A, j, i, val); alpar@1: } alpar@1: } alpar@1: fini: if (hbm != NULL) hbm_free_mat(hbm); alpar@1: return A; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_count_nnz - determine number of non-zeros in sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * int spm_count_nnz(const SPM *A); alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine spm_count_nnz returns the number of structural non-zero alpar@1: * elements in the specified sparse matrix. */ alpar@1: alpar@1: int spm_count_nnz(const SPM *A) alpar@1: { SPME *e; alpar@1: int i, nnz = 0; alpar@1: for (i = 1; i <= A->m; i++) alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) nnz++; alpar@1: return nnz; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_drop_zeros - remove zero elements from sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * int spm_drop_zeros(SPM *A, double eps); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_drop_zeros removes all elements from the specified alpar@1: * sparse matrix, whose absolute value is less than eps. alpar@1: * alpar@1: * If the parameter eps is 0, only zero elements are removed from the alpar@1: * matrix. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine returns the number of elements removed. */ alpar@1: alpar@1: int spm_drop_zeros(SPM *A, double eps) alpar@1: { SPME *e, *next; alpar@1: int i, count = 0; alpar@1: for (i = 1; i <= A->m; i++) alpar@1: { for (e = A->row[i]; e != NULL; e = next) alpar@1: { next = e->r_next; alpar@1: if (e->val == 0.0 || fabs(e->val) < eps) alpar@1: { /* remove element from the row list */ alpar@1: if (e->r_prev == NULL) alpar@1: A->row[e->i] = e->r_next; alpar@1: else alpar@1: e->r_prev->r_next = e->r_next; alpar@1: if (e->r_next == NULL) alpar@1: ; alpar@1: else alpar@1: e->r_next->r_prev = e->r_prev; alpar@1: /* remove element from the column list */ alpar@1: if (e->c_prev == NULL) alpar@1: A->col[e->j] = e->c_next; alpar@1: else alpar@1: e->c_prev->c_next = e->c_next; alpar@1: if (e->c_next == NULL) alpar@1: ; alpar@1: else alpar@1: e->c_next->c_prev = e->c_prev; alpar@1: /* return element to the memory pool */ alpar@1: dmp_free_atom(A->pool, e, sizeof(SPME)); alpar@1: count++; alpar@1: } alpar@1: } alpar@1: } alpar@1: return count; alpar@1: } alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_read_mat - read sparse matrix from text file alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_read_mat(const char *fname); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine reads a sparse matrix from a text file whose name is alpar@1: * specified by the parameter fname. alpar@1: * alpar@1: * For the file format see description of the routine spm_write_mat. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * On success the routine returns a pointer to the matrix created, alpar@1: * otherwise NULL. */ alpar@1: alpar@1: #if 1 alpar@1: SPM *spm_read_mat(const char *fname) alpar@1: { xassert(fname != fname); alpar@1: return NULL; alpar@1: } alpar@1: #else alpar@1: SPM *spm_read_mat(const char *fname) alpar@1: { SPM *A = NULL; alpar@1: PDS *pds; alpar@1: jmp_buf jump; alpar@1: int i, j, k, m, n, nnz, fail = 0; alpar@1: double val; alpar@1: xprintf("spm_read_mat: reading matrix from `%s'...\n", fname); alpar@1: pds = pds_open_file(fname); alpar@1: if (pds == NULL) alpar@1: { xprintf("spm_read_mat: unable to open `%s' - %s\n", fname, alpar@1: strerror(errno)); alpar@1: fail = 1; alpar@1: goto done; alpar@1: } alpar@1: if (setjmp(jump)) alpar@1: { fail = 1; alpar@1: goto done; alpar@1: } alpar@1: pds_set_jump(pds, jump); alpar@1: /* number of rows, number of columns, number of non-zeros */ alpar@1: m = pds_scan_int(pds); alpar@1: if (m < 0) alpar@1: pds_error(pds, "invalid number of rows\n"); alpar@1: n = pds_scan_int(pds); alpar@1: if (n < 0) alpar@1: pds_error(pds, "invalid number of columns\n"); alpar@1: nnz = pds_scan_int(pds); alpar@1: if (nnz < 0) alpar@1: pds_error(pds, "invalid number of non-zeros\n"); alpar@1: /* create matrix */ alpar@1: xprintf("spm_read_mat: %d rows, %d columns, %d non-zeros\n", alpar@1: m, n, nnz); alpar@1: A = spm_create_mat(m, n); alpar@1: /* read matrix elements */ alpar@1: for (k = 1; k <= nnz; k++) alpar@1: { /* row index, column index, element value */ alpar@1: i = pds_scan_int(pds); alpar@1: if (!(1 <= i && i <= m)) alpar@1: pds_error(pds, "row index out of range\n"); alpar@1: j = pds_scan_int(pds); alpar@1: if (!(1 <= j && j <= n)) alpar@1: pds_error(pds, "column index out of range\n"); alpar@1: val = pds_scan_num(pds); alpar@1: /* add new element to the matrix */ alpar@1: spm_new_elem(A, i, j, val); alpar@1: } alpar@1: xprintf("spm_read_mat: %d lines were read\n", pds->count); alpar@1: done: if (pds != NULL) pds_close_file(pds); alpar@1: if (fail && A != NULL) spm_delete_mat(A), A = NULL; alpar@1: return A; alpar@1: } alpar@1: #endif alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_write_mat - write sparse matrix to text file alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * int spm_write_mat(const SPM *A, const char *fname); alpar@1: * alpar@1: * DESCRIPTION alpar@1: * alpar@1: * The routine spm_write_mat writes the specified sparse matrix to a alpar@1: * text file whose name is specified by the parameter fname. This file alpar@1: * can be read back with the routine spm_read_mat. alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * On success the routine returns zero, otherwise non-zero. alpar@1: * alpar@1: * FILE FORMAT alpar@1: * alpar@1: * The file created by the routine spm_write_mat is a plain text file, alpar@1: * which contains the following information: alpar@1: * alpar@1: * m n nnz alpar@1: * row[1] col[1] val[1] alpar@1: * row[2] col[2] val[2] alpar@1: * . . . alpar@1: * row[nnz] col[nnz] val[nnz] alpar@1: * alpar@1: * where: alpar@1: * m is the number of rows; alpar@1: * n is the number of columns; alpar@1: * nnz is the number of non-zeros; alpar@1: * row[k], k = 1,...,nnz, are row indices; alpar@1: * col[k], k = 1,...,nnz, are column indices; alpar@1: * val[k], k = 1,...,nnz, are element values. */ alpar@1: alpar@1: #if 1 alpar@1: int spm_write_mat(const SPM *A, const char *fname) alpar@1: { xassert(A != A); alpar@1: xassert(fname != fname); alpar@1: return 0; alpar@1: } alpar@1: #else alpar@1: int spm_write_mat(const SPM *A, const char *fname) alpar@1: { FILE *fp; alpar@1: int i, nnz, ret = 0; alpar@1: xprintf("spm_write_mat: writing matrix to `%s'...\n", fname); alpar@1: fp = fopen(fname, "w"); alpar@1: if (fp == NULL) alpar@1: { xprintf("spm_write_mat: unable to create `%s' - %s\n", fname, alpar@1: strerror(errno)); alpar@1: ret = 1; alpar@1: goto done; alpar@1: } alpar@1: /* number of rows, number of columns, number of non-zeros */ alpar@1: nnz = spm_count_nnz(A); alpar@1: fprintf(fp, "%d %d %d\n", A->m, A->n, nnz); alpar@1: /* walk through rows of the matrix */ alpar@1: for (i = 1; i <= A->m; i++) alpar@1: { SPME *e; alpar@1: /* walk through elements of i-th row */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: { /* row index, column index, element value */ alpar@1: fprintf(fp, "%d %d %.*g\n", e->i, e->j, DBL_DIG, e->val); alpar@1: } alpar@1: } alpar@1: fflush(fp); alpar@1: if (ferror(fp)) alpar@1: { xprintf("spm_write_mat: writing error on `%s' - %s\n", fname, alpar@1: strerror(errno)); alpar@1: ret = 1; alpar@1: goto done; alpar@1: } alpar@1: xprintf("spm_write_mat: %d lines were written\n", 1 + nnz); alpar@1: done: if (fp != NULL) fclose(fp); alpar@1: return ret; alpar@1: } alpar@1: #endif alpar@1: alpar@1: /*********************************************************************** alpar@1: * NAME alpar@1: * alpar@1: * spm_transpose - transpose sparse matrix alpar@1: * alpar@1: * SYNOPSIS alpar@1: * alpar@1: * #include "glpspm.h" alpar@1: * SPM *spm_transpose(const SPM *A); alpar@1: * alpar@1: * RETURNS alpar@1: * alpar@1: * The routine computes and returns sparse matrix B, which is a matrix alpar@1: * transposed to sparse matrix A. */ alpar@1: alpar@1: SPM *spm_transpose(const SPM *A) alpar@1: { SPM *B; alpar@1: int i; alpar@1: B = spm_create_mat(A->n, A->m); alpar@1: for (i = 1; i <= A->m; i++) alpar@1: { SPME *e; alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: spm_new_elem(B, e->j, i, e->val); alpar@1: } alpar@1: return B; alpar@1: } alpar@1: alpar@1: SPM *spm_add_sym(const SPM *A, const SPM *B) alpar@1: { /* add two sparse matrices (symbolic phase) */ alpar@1: SPM *C; alpar@1: int i, j, *flag; alpar@1: xassert(A->m == B->m); alpar@1: xassert(A->n == B->n); alpar@1: /* create resultant matrix */ alpar@1: C = spm_create_mat(A->m, A->n); alpar@1: /* allocate and clear the flag array */ alpar@1: flag = xcalloc(1+C->n, sizeof(int)); alpar@1: for (j = 1; j <= C->n; j++) alpar@1: flag[j] = 0; alpar@1: /* compute pattern of C = A + B */ alpar@1: for (i = 1; i <= C->m; i++) alpar@1: { SPME *e; alpar@1: /* at the beginning i-th row of C is empty */ alpar@1: /* (i-th row of C) := (i-th row of C) union (i-th row of A) */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: { /* (note that i-th row of A may have duplicate elements) */ alpar@1: j = e->j; alpar@1: if (!flag[j]) alpar@1: { spm_new_elem(C, i, j, 0.0); alpar@1: flag[j] = 1; alpar@1: } alpar@1: } alpar@1: /* (i-th row of C) := (i-th row of C) union (i-th row of B) */ alpar@1: for (e = B->row[i]; e != NULL; e = e->r_next) alpar@1: { /* (note that i-th row of B may have duplicate elements) */ alpar@1: j = e->j; alpar@1: if (!flag[j]) alpar@1: { spm_new_elem(C, i, j, 0.0); alpar@1: flag[j] = 1; alpar@1: } alpar@1: } alpar@1: /* reset the flag array */ alpar@1: for (e = C->row[i]; e != NULL; e = e->r_next) alpar@1: flag[e->j] = 0; alpar@1: } alpar@1: /* check and deallocate the flag array */ alpar@1: for (j = 1; j <= C->n; j++) alpar@1: xassert(!flag[j]); alpar@1: xfree(flag); alpar@1: return C; alpar@1: } alpar@1: alpar@1: void spm_add_num(SPM *C, double alfa, const SPM *A, double beta, alpar@1: const SPM *B) alpar@1: { /* add two sparse matrices (numeric phase) */ alpar@1: int i, j; alpar@1: double *work; alpar@1: /* allocate and clear the working array */ alpar@1: work = xcalloc(1+C->n, sizeof(double)); alpar@1: for (j = 1; j <= C->n; j++) alpar@1: work[j] = 0.0; alpar@1: /* compute matrix C = alfa * A + beta * B */ alpar@1: for (i = 1; i <= C->n; i++) alpar@1: { SPME *e; alpar@1: /* work := alfa * (i-th row of A) + beta * (i-th row of B) */ alpar@1: /* (note that A and/or B may have duplicate elements) */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: work[e->j] += alfa * e->val; alpar@1: for (e = B->row[i]; e != NULL; e = e->r_next) alpar@1: work[e->j] += beta * e->val; alpar@1: /* (i-th row of C) := work, work := 0 */ alpar@1: for (e = C->row[i]; e != NULL; e = e->r_next) alpar@1: { j = e->j; alpar@1: e->val = work[j]; alpar@1: work[j] = 0.0; alpar@1: } alpar@1: } alpar@1: /* check and deallocate the working array */ alpar@1: for (j = 1; j <= C->n; j++) alpar@1: xassert(work[j] == 0.0); alpar@1: xfree(work); alpar@1: return; alpar@1: } alpar@1: alpar@1: SPM *spm_add_mat(double alfa, const SPM *A, double beta, const SPM *B) alpar@1: { /* add two sparse matrices (driver routine) */ alpar@1: SPM *C; alpar@1: C = spm_add_sym(A, B); alpar@1: spm_add_num(C, alfa, A, beta, B); alpar@1: return C; alpar@1: } alpar@1: alpar@1: SPM *spm_mul_sym(const SPM *A, const SPM *B) alpar@1: { /* multiply two sparse matrices (symbolic phase) */ alpar@1: int i, j, k, *flag; alpar@1: SPM *C; alpar@1: xassert(A->n == B->m); alpar@1: /* create resultant matrix */ alpar@1: C = spm_create_mat(A->m, B->n); alpar@1: /* allocate and clear the flag array */ alpar@1: flag = xcalloc(1+C->n, sizeof(int)); alpar@1: for (j = 1; j <= C->n; j++) alpar@1: flag[j] = 0; alpar@1: /* compute pattern of C = A * B */ alpar@1: for (i = 1; i <= C->m; i++) alpar@1: { SPME *e, *ee; alpar@1: /* compute pattern of i-th row of C */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: { k = e->j; alpar@1: for (ee = B->row[k]; ee != NULL; ee = ee->r_next) alpar@1: { j = ee->j; alpar@1: /* if a[i,k] != 0 and b[k,j] != 0 then c[i,j] != 0 */ alpar@1: if (!flag[j]) alpar@1: { /* c[i,j] does not exist, so create it */ alpar@1: spm_new_elem(C, i, j, 0.0); alpar@1: flag[j] = 1; alpar@1: } alpar@1: } alpar@1: } alpar@1: /* reset the flag array */ alpar@1: for (e = C->row[i]; e != NULL; e = e->r_next) alpar@1: flag[e->j] = 0; alpar@1: } alpar@1: /* check and deallocate the flag array */ alpar@1: for (j = 1; j <= C->n; j++) alpar@1: xassert(!flag[j]); alpar@1: xfree(flag); alpar@1: return C; alpar@1: } alpar@1: alpar@1: void spm_mul_num(SPM *C, const SPM *A, const SPM *B) alpar@1: { /* multiply two sparse matrices (numeric phase) */ alpar@1: int i, j; alpar@1: double *work; alpar@1: /* allocate and clear the working array */ alpar@1: work = xcalloc(1+A->n, sizeof(double)); alpar@1: for (j = 1; j <= A->n; j++) alpar@1: work[j] = 0.0; alpar@1: /* compute matrix C = A * B */ alpar@1: for (i = 1; i <= C->m; i++) alpar@1: { SPME *e, *ee; alpar@1: double temp; alpar@1: /* work := (i-th row of A) */ alpar@1: /* (note that A may have duplicate elements) */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: work[e->j] += e->val; alpar@1: /* compute i-th row of C */ alpar@1: for (e = C->row[i]; e != NULL; e = e->r_next) alpar@1: { j = e->j; alpar@1: /* c[i,j] := work * (j-th column of B) */ alpar@1: temp = 0.0; alpar@1: for (ee = B->col[j]; ee != NULL; ee = ee->c_next) alpar@1: temp += work[ee->i] * ee->val; alpar@1: e->val = temp; alpar@1: } alpar@1: /* reset the working array */ alpar@1: for (e = A->row[i]; e != NULL; e = e->r_next) alpar@1: work[e->j] = 0.0; alpar@1: } alpar@1: /* check and deallocate the working array */ alpar@1: for (j = 1; j <= A->n; j++) alpar@1: xassert(work[j] == 0.0); alpar@1: xfree(work); alpar@1: return; alpar@1: } alpar@1: alpar@1: SPM *spm_mul_mat(const SPM *A, const SPM *B) alpar@1: { /* multiply two sparse matrices (driver routine) */ alpar@1: SPM *C; alpar@1: C = spm_mul_sym(A, B); alpar@1: spm_mul_num(C, A, B); alpar@1: return C; alpar@1: } alpar@1: alpar@1: PER *spm_create_per(int n) alpar@1: { /* create permutation matrix */ alpar@1: PER *P; alpar@1: int k; alpar@1: xassert(n >= 0); alpar@1: P = xmalloc(sizeof(PER)); alpar@1: P->n = n; alpar@1: P->row = xcalloc(1+n, sizeof(int)); alpar@1: P->col = xcalloc(1+n, sizeof(int)); alpar@1: /* initially it is identity matrix */ alpar@1: for (k = 1; k <= n; k++) alpar@1: P->row[k] = P->col[k] = k; alpar@1: return P; alpar@1: } alpar@1: alpar@1: void spm_check_per(PER *P) alpar@1: { /* check permutation matrix for correctness */ alpar@1: int i, j; alpar@1: xassert(P->n >= 0); alpar@1: for (i = 1; i <= P->n; i++) alpar@1: { j = P->row[i]; alpar@1: xassert(1 <= j && j <= P->n); alpar@1: xassert(P->col[j] == i); alpar@1: } alpar@1: return; alpar@1: } alpar@1: alpar@1: void spm_delete_per(PER *P) alpar@1: { /* delete permutation matrix */ alpar@1: xfree(P->row); alpar@1: xfree(P->col); alpar@1: xfree(P); alpar@1: return; alpar@1: } alpar@1: alpar@1: /* eof */