SenK
SenK is a C++ library for high-performance linear solvers.
senk_io.hpp
Go to the documentation of this file.
1
7#ifndef SENK_IO_HPP
8#define SENK_IO_HPP
9
10#include <iostream>
11#include <fstream>
12#include <sstream>
13
14#include "senk_utils.hpp"
15#include "senk_helper.hpp"
16
17namespace senk {
18
22namespace io {
26enum Shape {
27 Sym,
28 Unsym
29};
42 std::string filename, double **val, int **cind, int **rptr,
43 int *N, int *M, Shape *shape, bool removeZeros)
44{
45 std::cout << "# Readfile MatrixMarket" << std::endl;
46 std::ifstream ifs(filename);
47 if(!ifs) {
48 std::cerr << "# Could not open input file." << std::endl;
49 return false;
50 }
51 std::string line;
52 std::getline(ifs, line);
53 if(line.find("complex") != std::string::npos) {
54 std::cerr << "# Header line is not valid." << std::endl;
55 return false;
56 }
57 if(line.find("symmetric") != std::string::npos) { shape[0] = Sym; }
58 else if(line.find("general") != std::string::npos) { shape[0] = Unsym; }
59 else { std::cerr << "# Header line is not valid." << std::endl; return false; }
60
61 int PE;
62 while(std::getline(ifs, line)) {
63 if(line[0] == '%') { continue; }
64 std::istringstream iss(line);
65 iss >> N[0] >> M[0] >> PE;
66 break;
67 }
68 std::cout << "# " << N[0] << " " << M[0] << " " << PE << std::endl;
69 *val = utils::SafeMalloc<double>(PE); //new double[PE];
70 *cind = utils::SafeMalloc<int>(PE); //new int[PE];
71 *rptr = utils::SafeMalloc<int>(N[0]+1); //new int[N[0]+1];
72 int *row = utils::SafeMalloc<int>(PE); //new int[PE];
73
74 int nnz = 0;
75 double t_val;
76 int t_row, t_col;
77 for(int i=0; i<PE; i++) {
78 ifs >> t_row >> t_col >> t_val;
79 if(removeZeros && t_val == 0) continue;
80 (*val)[nnz] = t_val;
81 row[nnz] = t_row-1;
82 (*cind)[nnz] = t_col-1;
83 nnz++;
84 }
85 if(nnz != PE) {
86 *val = utils::SafeRealloc<double>(*val, nnz);
87 *cind = utils::SafeRealloc<int>(*cind, nnz);
88 }
89 helper::QuickSort<int, int, double>(row, *cind, *val, 0, nnz-1);
90 int cnt = 0;
91 (*rptr)[0] = 0;
92 for(int i=0; i<N[0]; i++) {
93 int num = 0;
94 while(row[cnt] == i) {
95 num++; cnt++;
96 if(cnt == nnz) { break; }
97 }
98 (*rptr)[i+1] = (*rptr)[i] + num;
99 helper::QuickSort<int, double>(*cind, *val, (*rptr)[i], (*rptr)[i+1]-1);
100 }
101 std::free(row);
102 return true;
103}
104
105}
106
107}
108
109#endif
Shape
enum for the shape of matrices.
Definition: senk_io.hpp:26
bool ReadMatrixMarket(std::string filename, double **val, int **cind, int **rptr, int *N, int *M, Shape *shape, bool removeZeros)
Get a matrix in the CSR format from a MatrixMarket file.
Definition: senk_io.hpp:41
The top-level namespace of SenK.
Some supplemental functions are defined.
Utility functions are defined.