mmio.c 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010,2013-2017 CNRS
  4. * Copyright (C) 2008,2009,2011,2013,2016 Université de Bordeaux
  5. *
  6. * StarPU is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU Lesser General Public License as published by
  8. * the Free Software Foundation; either version 2.1 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * StarPU is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  14. *
  15. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  16. */
  17. /*
  18. * Matrix Market I/O library for ANSI C
  19. *
  20. * See http://math.nist.gov/MatrixMarket for details.
  21. *
  22. *
  23. */
  24. #include <stdio.h>
  25. #include <string.h>
  26. #include <stdlib.h>
  27. #include <ctype.h>
  28. #include "mmio.h"
  29. int mm_read_unsymmetric_sparse(const char *fname, int *M_, int *N_, int *nz_, double **val_, int **I_, int **J_)
  30. {
  31. FILE *f;
  32. MM_typecode matcode;
  33. int M, N, nz;
  34. int i;
  35. double *val;
  36. int *I, *J;
  37. if ((f = fopen(fname, "r")) == NULL)
  38. {
  39. fprintf(stderr, "File <%s> not found\n", fname);
  40. return -1;
  41. }
  42. if (mm_read_banner(f, &matcode) != 0)
  43. {
  44. fprintf(stderr, "mm_read_unsymetric: Could not process Matrix Market banner ");
  45. fprintf(stderr, " in file [%s]\n", fname);
  46. return -1;
  47. }
  48. if ( !(mm_is_real(matcode) && mm_is_matrix(matcode) && mm_is_sparse(matcode)))
  49. {
  50. fprintf(stderr, "Sorry, this application does not support ");
  51. fprintf(stderr, "Market Market type: [%s]\n", mm_typecode_to_str(matcode));
  52. return -1;
  53. }
  54. /* find out size of sparse matrix: M, N, nz .... */
  55. if (mm_read_mtx_crd_size(f, &M, &N, &nz) !=0)
  56. {
  57. fprintf(stderr, "read_unsymmetric_sparse(): could not parse matrix size.\n");
  58. return -1;
  59. }
  60. *M_ = M;
  61. *N_ = N;
  62. *nz_ = nz;
  63. /* reseve memory for matrices */
  64. I = (int *) malloc(nz * sizeof(int));
  65. J = (int *) malloc(nz * sizeof(int));
  66. val = (double *) malloc(nz * sizeof(double));
  67. *val_ = val;
  68. *I_ = I;
  69. *J_ = J;
  70. /* NOTE: when reading in doubles, ANSI C requires the use of the "l" */
  71. /* specifier as in "%lg", "%lf", "%le", otherwise errors will occur */
  72. /* (ANSI C X3.159-1989, Sec. 4.9.6.2, p. 136 lines 13-15) */
  73. for (i=0; i<nz; i++)
  74. {
  75. fscanf(f, "%d %d %lg\n", &I[i], &J[i], &val[i]);
  76. I[i]--; /* adjust from 1-based to 0-based */
  77. J[i]--;
  78. }
  79. fclose(f);
  80. return 0;
  81. }
  82. int mm_is_valid(MM_typecode matcode)
  83. {
  84. if (!mm_is_matrix(matcode)) return 0;
  85. if (mm_is_dense(matcode) && mm_is_pattern(matcode)) return 0;
  86. if (mm_is_real(matcode) && mm_is_hermitian(matcode)) return 0;
  87. if (mm_is_pattern(matcode) && (mm_is_hermitian(matcode) || mm_is_skew(matcode))) return 0;
  88. return 1;
  89. }
  90. int mm_read_banner(FILE *f, MM_typecode *matcode)
  91. {
  92. char line[MM_MAX_LINE_LENGTH+1];
  93. char banner[MM_MAX_TOKEN_LENGTH+1];
  94. char mtx[MM_MAX_TOKEN_LENGTH+1];
  95. char crd[MM_MAX_TOKEN_LENGTH+1];
  96. char data_type[MM_MAX_TOKEN_LENGTH+1];
  97. char storage_scheme[MM_MAX_TOKEN_LENGTH+1];
  98. char *p;
  99. mm_clear_typecode(matcode);
  100. if (fgets(line, MM_MAX_LINE_LENGTH, f) == NULL)
  101. return MM_PREMATURE_EOF;
  102. if (sscanf(line, "%"MM_MAX_TOKEN_LENGTH_S"s %"MM_MAX_TOKEN_LENGTH_S"s %"MM_MAX_TOKEN_LENGTH_S"s %"MM_MAX_TOKEN_LENGTH_S"s %"MM_MAX_TOKEN_LENGTH_S"s", banner, mtx, crd, data_type, storage_scheme) != 5)
  103. return MM_PREMATURE_EOF;
  104. for (p=mtx; *p!='\0'; *p=tolower(*p),p++); /* convert to lower case */
  105. for (p=crd; *p!='\0'; *p=tolower(*p),p++);
  106. for (p=data_type; *p!='\0'; *p=tolower(*p),p++);
  107. for (p=storage_scheme; *p!='\0'; *p=tolower(*p),p++);
  108. /* check for banner */
  109. if (strncmp(banner, MatrixMarketBanner, strlen(MatrixMarketBanner)) != 0)
  110. return MM_NO_HEADER;
  111. /* first field should be "mtx" */
  112. if (strcmp(mtx, MM_MTX_STR) != 0)
  113. return MM_UNSUPPORTED_TYPE;
  114. mm_set_matrix(matcode);
  115. /* second field describes whether this is a sparse matrix (in coordinate storage) or a dense array */
  116. if (strcmp(crd, MM_SPARSE_STR) == 0)
  117. mm_set_sparse(matcode);
  118. else if (strcmp(crd, MM_DENSE_STR) == 0)
  119. mm_set_dense(matcode);
  120. else
  121. return MM_UNSUPPORTED_TYPE;
  122. /* third field */
  123. if (strcmp(data_type, MM_REAL_STR) == 0)
  124. mm_set_real(matcode);
  125. else if (strcmp(data_type, MM_COMPLEX_STR) == 0)
  126. mm_set_complex(matcode);
  127. else if (strcmp(data_type, MM_PATTERN_STR) == 0)
  128. mm_set_pattern(matcode);
  129. else if (strcmp(data_type, MM_INT_STR) == 0)
  130. mm_set_integer(matcode);
  131. else
  132. return MM_UNSUPPORTED_TYPE;
  133. /* fourth field */
  134. if (strcmp(storage_scheme, MM_GENERAL_STR) == 0)
  135. mm_set_general(matcode);
  136. else if (strcmp(storage_scheme, MM_SYMM_STR) == 0)
  137. mm_set_symmetric(matcode);
  138. else if (strcmp(storage_scheme, MM_HERM_STR) == 0)
  139. mm_set_hermitian(matcode);
  140. else if (strcmp(storage_scheme, MM_SKEW_STR) == 0)
  141. mm_set_skew(matcode);
  142. else
  143. return MM_UNSUPPORTED_TYPE;
  144. return 0;
  145. }
  146. int mm_write_mtx_crd_size(FILE *f, int M, int N, int nz)
  147. {
  148. if (fprintf(f, "%d %d %d\n", M, N, nz) != 3)
  149. return MM_COULD_NOT_WRITE_FILE;
  150. else
  151. return 0;
  152. }
  153. int mm_read_mtx_crd_size(FILE *f, int *M, int *N, int *nz )
  154. {
  155. char line[MM_MAX_LINE_LENGTH];
  156. /* set return null parameter values, in case we exit with errors */
  157. *M = *N = *nz = 0;
  158. /* now continue scanning until you reach the end-of-comments */
  159. do
  160. {
  161. if (fgets(line,MM_MAX_LINE_LENGTH,f) == NULL)
  162. return MM_PREMATURE_EOF;
  163. } while (line[0] == '%');
  164. /* line[] is either blank or has M,N, nz */
  165. if (sscanf(line, "%d %d %d", M, N, nz) == 3)
  166. return 0;
  167. else
  168. {
  169. int num_items_read;
  170. do
  171. {
  172. num_items_read = fscanf(f, "%d %d %d", M, N, nz);
  173. if (num_items_read == EOF) return MM_PREMATURE_EOF;
  174. }
  175. while (num_items_read != 3);
  176. }
  177. return 0;
  178. }
  179. int mm_read_mtx_array_size(FILE *f, int *M, int *N)
  180. {
  181. char line[MM_MAX_LINE_LENGTH];
  182. /* set return null parameter values, in case we exit with errors */
  183. *M = *N = 0;
  184. /* now continue scanning until you reach the end-of-comments */
  185. do
  186. {
  187. if (fgets(line,MM_MAX_LINE_LENGTH,f) == NULL)
  188. return MM_PREMATURE_EOF;
  189. } while (line[0] == '%');
  190. /* line[] is either blank or has M,N, nz */
  191. if (sscanf(line, "%d %d", M, N) == 2)
  192. return 0;
  193. else /* we have a blank line */
  194. {
  195. int num_items_read;
  196. do
  197. {
  198. num_items_read = fscanf(f, "%d %d", M, N);
  199. if (num_items_read == EOF) return MM_PREMATURE_EOF;
  200. }
  201. while (num_items_read != 2);
  202. }
  203. return 0;
  204. }
  205. int mm_write_mtx_array_size(FILE *f, int M, int N)
  206. {
  207. if (fprintf(f, "%d %d\n", M, N) != 2)
  208. return MM_COULD_NOT_WRITE_FILE;
  209. else
  210. return 0;
  211. }
  212. /*-------------------------------------------------------------------------*/
  213. /******************************************************************/
  214. /* use when I[], J[], and val[]J, and val[] are already allocated */
  215. /******************************************************************/
  216. int mm_read_mtx_crd_data(FILE *f, int M, int N, int nz, int I[], int J[], double val[], MM_typecode matcode)
  217. {
  218. int i;
  219. if (mm_is_complex(matcode))
  220. {
  221. for (i=0; i<nz; i++)
  222. if (fscanf(f, "%d %d %lg %lg", &I[i], &J[i], &val[2*i], &val[2*i+1]) != 4)
  223. return MM_PREMATURE_EOF;
  224. }
  225. else if (mm_is_real(matcode))
  226. {
  227. for (i=0; i<nz; i++)
  228. {
  229. if (fscanf(f, "%d %d %lg\n", &I[i], &J[i], &val[i]) != 3)
  230. return MM_PREMATURE_EOF;
  231. }
  232. }
  233. else if (mm_is_pattern(matcode))
  234. {
  235. for (i=0; i<nz; i++)
  236. if (fscanf(f, "%d %d", &I[i], &J[i]) != 2)
  237. return MM_PREMATURE_EOF;
  238. }
  239. else
  240. return MM_UNSUPPORTED_TYPE;
  241. return 0;
  242. }
  243. int mm_read_mtx_crd_entry(FILE *f, int *I, int *J, double *real, double *imag, MM_typecode matcode)
  244. {
  245. if (mm_is_complex(matcode))
  246. {
  247. if (fscanf(f, "%d %d %lg %lg", I, J, real, imag) != 4)
  248. return MM_PREMATURE_EOF;
  249. }
  250. else if (mm_is_real(matcode))
  251. {
  252. if (fscanf(f, "%d %d %lg\n", I, J, real) != 3)
  253. return MM_PREMATURE_EOF;
  254. }
  255. else if (mm_is_pattern(matcode))
  256. {
  257. if (fscanf(f, "%d %d", I, J) != 2) return MM_PREMATURE_EOF;
  258. }
  259. else
  260. return MM_UNSUPPORTED_TYPE;
  261. return 0;
  262. }
  263. /************************************************************************
  264. mm_read_mtx_crd() fills M, N, nz, array of values, and return
  265. type code, e.g. 'MCRS'
  266. if matrix is complex, values[] is of size 2*nz,
  267. (nz pairs of real/imaginary values)
  268. ************************************************************************/
  269. int mm_read_mtx_crd(char *fname, int *M, int *N, int *nz, int **I, int **J, double **val, MM_typecode *matcode)
  270. {
  271. int ret_code;
  272. FILE *f;
  273. if (strcmp(fname, "stdin") == 0) f=stdin;
  274. else
  275. if ((f = fopen(fname, "r")) == NULL)
  276. return MM_COULD_NOT_READ_FILE;
  277. if ((ret_code = mm_read_banner(f, matcode)) != 0)
  278. {
  279. if (f != stdin) fclose(f);
  280. return ret_code;
  281. }
  282. if (!(mm_is_valid(*matcode) && mm_is_sparse(*matcode) && mm_is_matrix(*matcode)))
  283. {
  284. if (f != stdin) fclose(f);
  285. return MM_UNSUPPORTED_TYPE;
  286. }
  287. if ((ret_code = mm_read_mtx_crd_size(f, M, N, nz)) != 0)
  288. {
  289. if (f != stdin) fclose(f);
  290. return ret_code;
  291. }
  292. *I = (int *) malloc(*nz * sizeof(int));
  293. *J = (int *) malloc(*nz * sizeof(int));
  294. *val = NULL;
  295. if (mm_is_complex(*matcode))
  296. {
  297. *val = (double *) malloc(*nz * 2 * sizeof(double));
  298. ret_code = mm_read_mtx_crd_data(f, *M, *N, *nz, *I, *J, *val, *matcode);
  299. if (ret_code != 0)
  300. {
  301. if (f != stdin) fclose(f);
  302. return ret_code;
  303. }
  304. }
  305. else if (mm_is_real(*matcode))
  306. {
  307. *val = (double *) malloc(*nz * sizeof(double));
  308. ret_code = mm_read_mtx_crd_data(f, *M, *N, *nz, *I, *J, *val, *matcode);
  309. if (ret_code != 0)
  310. {
  311. if (f != stdin) fclose(f);
  312. return ret_code;
  313. }
  314. }
  315. else if (mm_is_pattern(*matcode))
  316. {
  317. ret_code = mm_read_mtx_crd_data(f, *M, *N, *nz, *I, *J, *val, *matcode);
  318. if (ret_code != 0)
  319. {
  320. if (f != stdin) fclose(f);
  321. return ret_code;
  322. }
  323. }
  324. if (f != stdin) fclose(f);
  325. return 0;
  326. }
  327. int mm_write_banner(FILE *f, MM_typecode matcode)
  328. {
  329. char *str = mm_typecode_to_str(matcode);
  330. int ret_code;
  331. ret_code = fprintf(f, "%s %s\n", MatrixMarketBanner, str);
  332. free(str);
  333. if (ret_code != 2)
  334. return MM_COULD_NOT_WRITE_FILE;
  335. else
  336. return 0;
  337. }
  338. int mm_write_mtx_crd(char fname[], int M, int N, int nz, int I[], int J[], double val[], MM_typecode matcode)
  339. {
  340. FILE *f;
  341. int i;
  342. if (strcmp(fname, "stdout") == 0)
  343. f = stdout;
  344. else if ((f = fopen(fname, "w")) == NULL)
  345. return MM_COULD_NOT_WRITE_FILE;
  346. /* print banner followed by typecode */
  347. fprintf(f, "%s ", MatrixMarketBanner);
  348. fprintf(f, "%s\n", mm_typecode_to_str(matcode));
  349. /* print matrix sizes and nonzeros */
  350. fprintf(f, "%d %d %d\n", M, N, nz);
  351. /* print values */
  352. if (mm_is_pattern(matcode))
  353. for (i=0; i<nz; i++)
  354. fprintf(f, "%d %d\n", I[i], J[i]);
  355. else if (mm_is_real(matcode))
  356. for (i=0; i<nz; i++)
  357. fprintf(f, "%d %d %20.16g\n", I[i], J[i], val[i]);
  358. else if (mm_is_complex(matcode))
  359. for (i=0; i<nz; i++)
  360. fprintf(f, "%d %d %20.16g %20.16g\n", I[i], J[i], val[2*i], val[2*i+1]);
  361. else
  362. {
  363. if (f != stdout) fclose(f);
  364. return MM_UNSUPPORTED_TYPE;
  365. }
  366. if (f !=stdout) fclose(f);
  367. return 0;
  368. }
  369. /**
  370. * Create a new copy of a string s. mm_strdup() is a common routine, but
  371. * not part of ANSI C, so it is included here. Used by mm_typecode_to_str().
  372. *
  373. */
  374. char *mm_strdup(const char *s)
  375. {
  376. int len = strlen(s);
  377. char *s2 = (char *) malloc((len+1)*sizeof(char));
  378. return strcpy(s2, s);
  379. }
  380. char *mm_typecode_to_str(MM_typecode matcode)
  381. {
  382. char buffer[MM_MAX_LINE_LENGTH];
  383. char *types[4];
  384. /* char *mm_strdup(const char *); */
  385. /* check for MTX type */
  386. if (mm_is_matrix(matcode))
  387. types[0] = MM_MTX_STR;
  388. /* check for CRD or ARR matrix */
  389. if (mm_is_sparse(matcode))
  390. types[1] = MM_SPARSE_STR;
  391. else if (mm_is_dense(matcode))
  392. types[1] = MM_DENSE_STR;
  393. else
  394. return NULL;
  395. /* check for element data type */
  396. if (mm_is_real(matcode))
  397. types[2] = MM_REAL_STR;
  398. else if (mm_is_complex(matcode))
  399. types[2] = MM_COMPLEX_STR;
  400. else if (mm_is_pattern(matcode))
  401. types[2] = MM_PATTERN_STR;
  402. else if (mm_is_integer(matcode))
  403. types[2] = MM_INT_STR;
  404. else
  405. return NULL;
  406. /* check for symmetry type */
  407. if (mm_is_general(matcode))
  408. types[3] = MM_GENERAL_STR;
  409. else if (mm_is_symmetric(matcode))
  410. types[3] = MM_SYMM_STR;
  411. else if (mm_is_hermitian(matcode))
  412. types[3] = MM_HERM_STR;
  413. else if (mm_is_skew(matcode))
  414. types[3] = MM_SKEW_STR;
  415. else
  416. return NULL;
  417. snprintf(buffer, sizeof(buffer), "%s %s %s %s", types[0], types[1], types[2], types[3]);
  418. return mm_strdup(buffer);
  419. }