SpM Handbook 1.2.4
All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
spm_io.c
Go to the documentation of this file.
1/**
2 *
3 * @file spm_io.c
4 *
5 * SParse Matrix package I/O routines.
6 *
7 * @copyright 2016-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8 * Univ. Bordeaux. All rights reserved.
9 *
10 * @version 1.2.4
11 * @author Pierre Ramet
12 * @author Mathieu Faverge
13 * @author Matias Hastaran
14 * @author Tony Delarue
15 * @author Alycia Lisito
16 * @date 2024-06-26
17 *
18 **/
19#include "common.h"
20
21/**
22 *******************************************************************************
23 *
24 * @ingroup spm_dev_io
25 *
26 * @brief Read an array of integer.
27 *
28 *******************************************************************************
29 *
30 * @param[in] stream
31 * The opened file in which the spm is stored.
32 *
33 * @param[in] n
34 * Number of elements.
35 *
36 * @param[out] array
37 * Array of results.
38 *
39 ********************************************************************************
40 *
41 * @retval SPM_SUCCESS if the read happened successfully,
42 * @retval SPM_ERR_FILE if the input format is incorrect.
43 *
44 *******************************************************************************/
45static inline int
46readArrayOfInteger( FILE *stream,
47 spm_int_t n,
48 spm_int_t *array )
49{
50 long tmp1, tmp2, tmp3, tmp4;
51 spm_int_t i;
52
53 /* Read 4 by 4 */
54 for (i=0; i<(n-3); i+=4)
55 {
56 if (4 != fscanf(stream, "%ld %ld %ld %ld", &tmp1, &tmp2, &tmp3, &tmp4)){
57 spm_print_error("spmLoad: Wrong input format");
58 return SPM_ERR_FILE;
59 }
60
61 array[i ] = (spm_int_t)tmp1;
62 array[i+1] = (spm_int_t)tmp2;
63 array[i+2] = (spm_int_t)tmp3;
64 array[i+3] = (spm_int_t)tmp4;
65 }
66
67 assert( n-i < 4 );
68 switch ( n - i )
69 {
70 case 3:
71 if (3 != fscanf(stream, "%ld %ld %ld", &tmp1, &tmp2, &tmp3)){
72 spm_print_error("spmLoad: Wrong input format");
73 return SPM_ERR_FILE;
74 }
75
76 array[i ] = (spm_int_t)tmp1;
77 array[i+1] = (spm_int_t)tmp2;
78 array[i+2] = (spm_int_t)tmp3;
79 break;
80 case 2:
81 if (2 != fscanf(stream, "%ld %ld", &tmp1, &tmp2)){
82 spm_print_error("spmLoad: Wrong input format");
83 return SPM_ERR_FILE;
84 }
85
86 array[i ] = (spm_int_t)tmp1;
87 array[i+1] = (spm_int_t)tmp2;
88 break;
89 case 1:
90 if (1 != fscanf(stream, "%ld", &tmp1)){
91 spm_print_error("spmLoad: Wrong input format");
92 return SPM_ERR_FILE;
93 }
94
95 array[i ] = (spm_int_t)tmp1;
96 break;
97 }
98
99 return SPM_SUCCESS;
100}
101
102/**
103 *******************************************************************************
104 *
105 * @ingroup spm_dev_io
106 *
107 * @brief Read an array of 64bits complex.
108 *
109 *******************************************************************************
110 *
111 * @param[in] stream
112 * The opened file in which the spm is stored.
113 *
114 * @param[in] n
115 * Number of elements.
116 *
117 * @param[out] array
118 * Array of results.
119 *
120 ********************************************************************************
121 *
122 * @retval SPM_SUCCESS if the read happened successfully,
123 * @retval SPM_ERR_FILE if the input format is incorrect.
124 *
125 *******************************************************************************/
126static inline int
128 spm_int_t n,
129 spm_complex64_t *array )
130{
131 double tmp1, tmp2, tmp3, tmp4;
132 double tmp5, tmp6, tmp7, tmp8;
133 spm_int_t i;
134
135 /* Read 4 by 4 */
136 for (i=0; i<(n-3); i+=4)
137 {
138 if (8 != fscanf(stream, "%lg %lg %lg %lg %lg %lg %lg %lg",
139 &tmp1, &tmp2, &tmp3, &tmp4,
140 &tmp5, &tmp6, &tmp7, &tmp8)){
141 spm_print_error("spmLoad: Wrong input format");
142 return SPM_ERR_FILE;
143 }
144 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
145 array[i+1] = (spm_complex64_t)(tmp3 + I * tmp4);
146 array[i+2] = (spm_complex64_t)(tmp5 + I * tmp6);
147 array[i+3] = (spm_complex64_t)(tmp7 + I * tmp8);
148 }
149
150 assert( n-i < 4 );
151 switch ( n - i )
152 {
153 case 3:
154 if (6 != fscanf(stream, "%lg %lg %lg %lg %lg %lg",
155 &tmp1, &tmp2, &tmp3, &tmp4,
156 &tmp5, &tmp6)){
157 spm_print_error("spmLoad: Wrong input format");
158 return SPM_ERR_FILE;
159 }
160 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
161 array[i+1] = (spm_complex64_t)(tmp3 + I * tmp4);
162 array[i+2] = (spm_complex64_t)(tmp5 + I * tmp6);
163 break;
164
165 case 2:
166 if (4 != fscanf(stream, "%lg %lg %lg %lg",
167 &tmp1, &tmp2, &tmp3, &tmp4)){
168 spm_print_error("spmLoad: Wrong input format");
169 return SPM_ERR_FILE;
170 }
171 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
172 array[i+1] = (spm_complex64_t)(tmp3 + I * tmp4);
173 break;
174
175 case 1:
176 if (2 != fscanf(stream, "%lg %lg",
177 &tmp1, &tmp2)){
178 spm_print_error("spmLoad: Wrong input format");
179 return SPM_ERR_FILE;
180 }
181 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
182 break;
183 }
184
185 return SPM_SUCCESS;
186}
187
188/**
189 *******************************************************************************
190 *
191 * @ingroup spm_dev_io
192 *
193 * @brief Read an array of 32bits complex.
194 *
195 *******************************************************************************
196 *
197 * @param[in] stream
198 * The opened file in which the spm is stored.
199 *
200 * @param[in] n
201 * Number of elements.
202 *
203 * @param[out] array
204 * Array of results.
205 *
206 ********************************************************************************
207 *
208 * @retval SPM_SUCCESS if the read happened successfully,
209 * @retval SPM_ERR_FILE if the input format is incorrect.
210 *
211 *******************************************************************************/
212static inline int
214 spm_int_t n,
215 spm_complex32_t *array )
216{
217 float tmp1, tmp2, tmp3, tmp4;
218 float tmp5, tmp6, tmp7, tmp8;
219 spm_int_t i;
220
221 /* Read 4 by 4 */
222 for (i=0; i<(n-3); i+=4)
223 {
224 if (8 != fscanf(stream, "%g %g %g %g %g %g %g %g",
225 &tmp1, &tmp2, &tmp3, &tmp4,
226 &tmp5, &tmp6, &tmp7, &tmp8)){
227 spm_print_error("spmLoad: Wrong input format");
228 return SPM_ERR_FILE;
229 }
230 array[i ] = (spm_complex32_t)(tmp1 + I * tmp2);
231 array[i+1] = (spm_complex32_t)(tmp3 + I * tmp4);
232 array[i+2] = (spm_complex32_t)(tmp5 + I * tmp6);
233 array[i+3] = (spm_complex32_t)(tmp7 + I * tmp8);
234 }
235
236 assert( n-i < 4 );
237 switch ( n - i )
238 {
239 case 3:
240 if (6 != fscanf(stream, "%g %g %g %g %g %g",
241 &tmp1, &tmp2, &tmp3, &tmp4,
242 &tmp5, &tmp6)){
243 spm_print_error("spmLoad: Wrong input format");
244 return SPM_ERR_FILE;
245 }
246 array[i ] = (spm_complex32_t)(tmp1 + I * tmp2);
247 array[i+1] = (spm_complex32_t)(tmp3 + I * tmp4);
248 array[i+2] = (spm_complex32_t)(tmp5 + I * tmp6);
249 break;
250
251 case 2:
252 if (4 != fscanf(stream, "%g %g %g %g",
253 &tmp1, &tmp2, &tmp3, &tmp4)){
254 spm_print_error("spmLoad: Wrong input format");
255 return SPM_ERR_FILE;
256 }
257 array[i ] = (spm_complex32_t)(tmp1 + I * tmp2);
258 array[i+1] = (spm_complex32_t)(tmp3 + I * tmp4);
259 break;
260
261 case 1:
262 if (2 != fscanf(stream, "%g %g",
263 &tmp1, &tmp2)){
264 spm_print_error("spmLoad: Wrong input format");
265 return SPM_ERR_FILE;
266 }
267 array[i ] = (spm_complex32_t)(tmp1 + I * tmp2);
268 break;
269 }
270
271 return SPM_SUCCESS;
272}
273
274/**
275 *******************************************************************************
276 *
277 * @ingroup spm_dev_io
278 *
279 * @brief Read an array of double.
280 *
281 *******************************************************************************
282 *
283 * @param[in] stream
284 * The opened file in which the spm is stored.
285 *
286 * @param[in] n
287 * Number of elements.
288 *
289 * @param[out] array
290 * Array of results.
291 *
292 ********************************************************************************
293 *
294 * @retval SPM_SUCCESS if the read happened successfully,
295 * @retval SPM_ERR_FILE if the input format is incorrect.
296 *
297 *******************************************************************************/
298static inline int
299readArrayOfDouble( FILE *stream,
300 spm_int_t n,
301 double *array )
302{
303 double tmp1, tmp2, tmp3, tmp4;
304 spm_int_t i;
305
306 /* Read 4 by 4 */
307 for (i=0; i<(n-3); i+=4)
308 {
309 if (4 != fscanf(stream, "%lg %lg %lg %lg",
310 &tmp1, &tmp2, &tmp3, &tmp4)){
311 spm_print_error("spmLoad: Wrong input format");
312 return SPM_ERR_FILE;
313 }
314 array[i ] = (double)(tmp1);
315 array[i+1] = (double)(tmp2);
316 array[i+2] = (double)(tmp3);
317 array[i+3] = (double)(tmp4);
318 }
319
320 assert( n-i < 4 );
321 switch ( n - i )
322 {
323 case 3:
324 if (1 != fscanf(stream, "%lg %lg %lg",
325 &tmp1, &tmp2, &tmp3)){
326 spm_print_error("spmLoad: Wrong input format");
327 return SPM_ERR_FILE;
328 }
329 array[i ] = (double)(tmp1);
330 array[i+1] = (double)(tmp2);
331 array[i+2] = (double)(tmp3);
332 break;
333
334 case 2:
335 if (2 != fscanf(stream, "%lg %lg",
336 &tmp1, &tmp2)){
337 spm_print_error("spmLoad: Wrong input format");
338 return SPM_ERR_FILE;
339 }
340 array[i ] = (double)(tmp1);
341 array[i+1] = (double)(tmp2);
342 break;
343
344 case 1:
345 if (1 != fscanf(stream, "%lg",
346 &tmp1)){
347 spm_print_error("spmLoad: Wrong input format");
348 return SPM_ERR_FILE;
349 }
350 array[i ] = (double)(tmp1);
351 break;
352 }
353
354 return SPM_SUCCESS;
355}
356
357
358/**
359 *******************************************************************************
360 *
361 * @ingroup spm_dev_io
362 *
363 * @brief Read an array of float.
364 *
365 *******************************************************************************
366 *
367 * @param[in] stream
368 * The opened file in which the spm is stored.
369 *
370 * @param[in] n
371 * Number of elements.
372 *
373 * @param[out] array
374 * Array of results.
375 *
376 ********************************************************************************
377 *
378 * @retval SPM_SUCCESS if the read happened successfully,
379 * @retval SPM_ERR_FILE if the input format is incorrect.
380 *
381 *******************************************************************************/
382static inline int
383readArrayOfFloat( FILE *stream,
384 spm_int_t n,
385 float *array )
386{
387 float tmp1, tmp2, tmp3, tmp4;
388 spm_int_t i;
389
390 /* Read 4 by 4 */
391 for (i=0; i<(n-3); i+=4)
392 {
393 if (4 != fscanf(stream, "%g %g %g %g",
394 &tmp1, &tmp2, &tmp3, &tmp4)){
395 spm_print_error("spmLoad: Wrong input format");
396 return SPM_ERR_FILE;
397 }
398 array[i ] = (float)(tmp1);
399 array[i+1] = (float)(tmp2);
400 array[i+2] = (float)(tmp3);
401 array[i+3] = (float)(tmp4);
402 }
403
404 assert( n-i < 4 );
405 switch ( n - i )
406 {
407 case 3:
408 if (3 != fscanf(stream, "%g %g %g",
409 &tmp1, &tmp2, &tmp3)){
410 spm_print_error("spmLoad: Wrong input format");
411 return SPM_ERR_FILE;
412 }
413 array[i ] = (float)(tmp1);
414 array[i+1] = (float)(tmp2);
415 array[i+2] = (float)(tmp3);
416 break;
417
418 case 2:
419 if (2 != fscanf(stream, "%g %g",
420 &tmp1, &tmp2)){
421 spm_print_error("spmLoad: Wrong input format");
422 return SPM_ERR_FILE;
423 }
424 array[i ] = (float)(tmp1);
425 array[i+1] = (float)(tmp2);
426 break;
427
428 case 1:
429 if (1 != fscanf(stream, "%g", &tmp1)){
430 spm_print_error("spmLoad: Wrong input format");
431 return SPM_ERR_FILE;
432 }
433 array[i ] = (float)(tmp1);
434 break;
435 }
436
437 return SPM_SUCCESS;
438}
439
440/**
441 *******************************************************************************
442 *
443 * @ingroup spm
444 *
445 * @brief Load the spm structure from a file (internal format). One node only.
446 *
447 *******************************************************************************
448 *
449 * @param[inout] spm
450 * On entry, an allocated spm data structure.
451 * On exit, the spm filled with the information read in the file.
452 *
453 * @param[in] filename
454 * The opened file in which the spm is stored. If filename == NULL,
455 * matrix.spm is opened.
456 *
457 *******************************************************************************
458 *
459 * @retval SPM_SUCCESS if the spm was load successfully.
460 * @retval SPM_ERR_FILE if the input format is incorrect.
461 *
462 *******************************************************************************/
463static inline int
465 const char *filename )
466{
467 FILE *infile = NULL;
468 spm_int_t colsize = 0;
469 spm_int_t rowsize = 0;
470 int rc = SPM_SUCCESS;
471 char line[256];
472 char *test;
473 long gN, nnzexp;
474
475 if ( filename == NULL ) {
476 filename = "matrix.spm";
477 }
478
479 infile = fopen( filename, "r" );
480 if ( infile == NULL ) {
481 spm_print_error( "spmLoad: Impossible to open the file matrix.spm\n");
482 return SPM_ERR_FILE;
483 }
484
485 /*
486 * Skip comments
487 */
488 do {
489 test = fgets( line, 256, infile );
490 if ( test != line ) {
491 fclose( infile );
492 return SPM_ERR_FILE;
493 }
494 }
495 while( line[0] == '#' );
496
497 /* Cleanup the line for coverity */
498 {
499 int i;
500 for (i=0; i<256; i++) {
501 if ( (line[i] == EOF ) || (line[i] == '\n') )
502 {
503 line[i] = '\0';
504 }
505 }
506 }
507
508 /*
509 * Read header
510 */
511 {
512 int version, mtxtype, flttype, fmttype, dof, layout;
513 long n, nnz;
514
515 if ( 10 != sscanf( line, "%d %d %d %d %ld %ld %ld %d %ld %d\n",
516 &version, &mtxtype, &flttype, &fmttype,
517 &gN, &n, &nnz, &dof, &nnzexp, &layout ) )
518 {
519 fclose( infile );
520 return SPM_ERR_FILE;
521 }
522
523 if ( ( gN <= 0 ) || ( n <= 0 ) || ( nnz <= 0 ) )
524 {
525 fclose( infile );
526 return SPM_ERR_FILE;
527 }
528
529 /* Handle only version 1 for now */
530 if (version != 1) {
531 fclose( infile );
532 return SPM_ERR_FILE;
533 }
534
535 spm->replicated = 1;
536 spm->mtxtype = (spm_mtxtype_t)mtxtype;
537 spm->flttype = (spm_coeftype_t)flttype;
538 spm->fmttype = (spm_fmttype_t)fmttype;
539 spm->gN = gN;
540 spm->n = n;
541 spm->nnz = nnz;
542 spm->dof = dof;
543 spm->layout = (spm_layout_t)layout;
544 }
545
546 switch(spm->fmttype){
547 case SpmCSC:
548 colsize = spm->n + 1;
549 rowsize = spm->nnz;
550 break;
551 case SpmCSR:
552 colsize = spm->nnz;
553 rowsize = spm->n + 1;
554 break;
555 case SpmIJV:
556 colsize = spm->nnz;
557 rowsize = spm->nnz;
558 break;
559 }
560
561 /*
562 * Read colptr
563 */
564 spm->colptr = malloc( colsize * sizeof(spm_int_t) );
565 rc = readArrayOfInteger( infile, colsize, spm->colptr );
566 if (rc != SPM_SUCCESS ) {
567 fclose( infile );
568 spmExit( spm );
569 return rc;
570 }
571
572 /*
573 * Read rowptr
574 */
575 spm->rowptr = malloc( rowsize * sizeof(spm_int_t) );
576 rc = readArrayOfInteger( infile, rowsize, spm->rowptr );
577 if (rc != SPM_SUCCESS ) {
578 fclose( infile );
579 spmExit( spm );
580 return rc;
581 }
582
583 /*
584 * Read dofs
585 */
586 if ( spm->dof > 0 ) {
587 spm->dofs = NULL;
588 }
589 else {
590 spm->dofs = malloc( (spm->n+1) * sizeof(spm_int_t) );
591 rc = readArrayOfInteger( infile, spm->n+1, spm->dofs );
592 if (rc != SPM_SUCCESS ) {
593 fclose( infile );
594 spmExit( spm );
595 return rc;
596 }
597 }
598
599 /* Update fields after reading the dof array to be variadic compatible */
601
602 /* Check that the computed values match the ones read */
603 assert( nnzexp == spm->nnzexp );
604 assert( spm->gN == gN );
605
606 /*
607 * Read values
608 */
609 if (spm->flttype == SpmPattern ) {
610 spm->values = NULL;
611 }
612 else {
613 spm->values = malloc( spm->nnzexp * spm_size_of( spm->flttype ) );
614 }
615
616 switch( spm->flttype ) {
617 case SpmPattern:
618 break;
619 case SpmFloat:
620 rc = readArrayOfFloat( infile, spm->nnzexp, spm->values );
621 break;
622 case SpmDouble:
623 rc = readArrayOfDouble( infile, spm->nnzexp, spm->values );
624 break;
625 case SpmComplex32:
626 rc = readArrayOfComplex32( infile, spm->nnzexp, spm->values );
627 break;
628 case SpmComplex64:
629 rc = readArrayOfComplex64( infile, spm->nnzexp, spm->values );
630 break;
631 }
632
633 fclose(infile);
634 return SPM_SUCCESS;
635}
636
637/**
638 *******************************************************************************
639 *
640 * @ingroup spm
641 *
642 * @brief Load the spm structure from a file (internal format).
643 *
644 *******************************************************************************
645 *
646 * @param[inout] spm
647 * On entry, an allocated spm data structure.
648 * On exit, the spm filled with the information read in the file.
649 *
650 * @param[in] filename
651 * The filename where the spm is stored. If filename == NULL,
652 * "./matrix.spm" is used.
653 *
654 * @param[in] comm
655 * The MPI communicator which share the file to load.
656 *
657 *******************************************************************************
658 *
659 * @retval SPM_SUCCESS if the load happened successfully,
660 * @retval SPM_ERR_FILE if the input format is incorrect.
661 *
662 *******************************************************************************/
663int
665 const char *filename,
666 SPM_Comm comm )
667{
668 spmatrix_t spmlocal;
669 int rc = SPM_SUCCESS;
670 int clustnum = 0;
671 spm_only_with_mpi int clustnbr = 1;
672 spm_only_with_mpi int distbycol = -1;
673
674#if defined(SPM_WITH_MPI)
675 /* Init to get the rank */
676 MPI_Comm_rank( comm, &clustnum );
677 MPI_Comm_size( comm, &clustnbr );
678#endif
679
680 /* Init the spm to know the rank */
681 if ( clustnum == 0 ) {
682 spmInitDist( &spmlocal, MPI_COMM_SELF );
683 rc = spm_load_local( &spmlocal, filename );
684 distbycol = (spmlocal.fmttype == SpmCSR) ? 0 : 1;
685 }
686
687#if defined(SPM_WITH_MPI)
688 MPI_Bcast( &rc, 1, MPI_INT, 0, comm );
689#endif
690
691 if ( rc != SPM_SUCCESS ) {
692 return rc;
693 }
694
695#if defined(SPM_WITH_MPI)
696 /* Scatter the spm if multiple nodes */
697 if ( clustnbr > 1 )
698 {
699 spmatrix_t spmdist;
700
701 /* Scatter the spm to all the process */
702 spmScatter( &spmdist, 0,
703 (clustnum == 0) ? &spmlocal : NULL,
704 0, NULL,
705 distbycol,
706 comm );
707
708 /* Switch the spm */
709 if ( clustnum == 0 ) {
710 spmExit( &spmlocal );
711 }
712 memcpy( spm, &spmdist, sizeof(spmatrix_t) );
713 }
714 else
715#endif
716 {
717 assert( clustnum == 0 );
718 memcpy( spm, &spmlocal, sizeof(spmatrix_t) );
719 }
720
721 (void)comm;
722 return SPM_SUCCESS;
723}
724
725/**
726 *******************************************************************************
727 *
728 * @ingroup spm
729 *
730 * @brief Load the spm structure from a file (internal format).
731 *
732 *******************************************************************************
733 *
734 * @param[inout] spm
735 * On entry, an allocated spm data structure.
736 * On exit, the spm filled with the information read in the file.
737 *
738 * @param[in] filename
739 * The filename where the spm is stored. If filename == NULL,
740 * "./matrix.spm" is used.
741 *
742 *******************************************************************************
743 *
744 * @retval SPM_SUCCESS if the load happened successfully,
745 * @retval SPM_ERR_FILE if the input format is incorrect.
746 *
747 *******************************************************************************/
748int
750 const char *filename )
751{
752 return spmLoadDist( spm, filename, MPI_COMM_WORLD );
753}
754
755/**
756 *******************************************************************************
757 *
758 * @ingroup spm_dev_io
759 *
760 * @brief write an array of 64bits complex.
761 *
762 *******************************************************************************
763 *
764 * @param[in] outfile
765 * The opened file in which to store the spm.
766 *
767 * @param[in] n
768 * numbers of elements.
769 *
770 * @param[in] array
771 * array to write.
772 *
773 *******************************************************************************
774 *
775 * @return SPM_SUCCESS if the write happened successfully.
776 *
777 *******************************************************************************/
778static inline int
780 spm_int_t n,
781 const spm_complex64_t *array )
782{
783 spm_int_t i;
784
785 /* Write 4 by 4 */
786 for (i=0; i<n; i++)
787 {
788 fprintf(outfile, "%e %e ", creal(array[i]), cimag(array[i]));
789 if (i%4 == 3) fprintf(outfile, "\n");
790 }
791 if ((i-1)%4 !=3) fprintf(outfile, "\n");
792
793 return SPM_SUCCESS;
794}
795
796/**
797 *******************************************************************************
798 *
799 * @ingroup spm_dev_io
800 *
801 * @brief write an array of 32bits complex.
802 *
803 *******************************************************************************
804 *
805 * @param[in] outfile
806 * The opened file in which to store the spm.
807 *
808 * @param[in] n
809 * numbers of elements.
810 *
811 * @param[in] array
812 * array to write.
813 *
814 *******************************************************************************
815 *
816 * @return SPM_SUCCESS if the write happened successfully.
817 *
818 *******************************************************************************/
819static inline int
821 spm_int_t n,
822 const spm_complex32_t *array )
823{
824 spm_int_t i;
825
826 /* Write 4 by 4 */
827 for (i=0; i<n; i++)
828 {
829 fprintf(outfile, "%e %e ", crealf(array[i]), cimagf(array[i]));
830 if (i%4 == 3) fprintf(outfile, "\n");
831 }
832 if ((i-1)%4 !=3) fprintf(outfile, "\n");
833
834 return SPM_SUCCESS;
835}
836
837/**
838 *******************************************************************************
839 *
840 * @ingroup spm_dev_io
841 *
842 * @brief write an array of double.
843 *
844 *******************************************************************************
845 *
846 * @param[in] outfile
847 * The opened file in which to store the spm.
848 *
849 * @param[in] n
850 * numbers of elements.
851 *
852 * @param[in] array
853 * array to write.
854 *
855 *******************************************************************************
856 *
857 * @return SPM_SUCCESS if the write happened successfully.
858 *
859 *******************************************************************************/
860static inline int
861writeArrayOfDouble( FILE *outfile,
862 spm_int_t n,
863 const double *array )
864{
865 spm_int_t i;
866
867 /* Write 4 by 4 */
868 for (i=0; i<n; i++)
869 {
870 fprintf(outfile, "%e ", array[i]);
871 if (i%4 == 3) fprintf(outfile, "\n");
872 }
873 if ((i-1)%4 !=3) fprintf(outfile, "\n");
874
875 return SPM_SUCCESS;
876}
877
878/**
879 *******************************************************************************
880 *
881 * @ingroup spm_dev_io
882 *
883 * @brief write an array of float.
884 *
885 *******************************************************************************
886 *
887 * @param[in] outfile
888 * The opened file in which to store the spm.
889 *
890 * @param[in] n
891 * numbers of elements.
892 *
893 * @param[in] array
894 * array to write.
895 *
896 *******************************************************************************
897 *
898 * @return SPM_SUCCESS if the write happened successfully.
899 *
900 *******************************************************************************/
901static inline int
902writeArrayOfFloat( FILE *outfile,
903 spm_int_t n,
904 const float *array )
905{
906 spm_int_t i;
907
908 /* Write 4 by 4 */
909 for (i=0; i<n; i++)
910 {
911 fprintf(outfile, "%e ", array[i]);
912 if (i%4 == 3) fprintf(outfile, "\n");
913 }
914 if ((i-1)%4 !=3) fprintf(outfile, "\n");
915
916 return SPM_SUCCESS;
917}
918
919/**
920 *******************************************************************************
921 *
922 * @ingroup spm
923 *
924 * @brief Save the global spm structure into a file (internal format).
925 *
926 *******************************************************************************
927 *
928 * @param[in] spm
929 * The sparse matrix to write into the file.
930 *
931 * @param[in] filename
932 * The filename in which to store the spm. If filename == NULL, data
933 * is saved into matrix.spm file.
934 *
935 ********************************************************************************
936 *
937 * @retval SPM_SUCCESS if the save happened successfully.
938 * @retval SPM_ERR_FILE if the input format is incorrect.
939 *
940 *******************************************************************************/
941static inline int
943 const char *filename )
944{
945 FILE *outfile = NULL;
946 spm_int_t i, colsize, rowsize;
947
948 if ( filename == NULL ) {
949 filename = "matrix.spm";
950 }
951
952 outfile = fopen( filename, "w" );
953 if ( outfile == NULL ) {
954 spm_print_error( "spmSave: Impossible to open the file matrix.spm\n");
955 return SPM_ERR_FILE;
956 }
957
958 /*
959 * Write header
960 */
961 fprintf( outfile,
962 "# version mtxtype flttype fmttype gN n nnz dof nnzexp layout\n"
963 "%d %d %d %d %ld %ld %ld %d %ld %d\n",
964 1, spm->mtxtype, spm->flttype, spm->fmttype,
965 (long)spm->gN, (long)spm->n, (long)spm->nnz,
966 (int)spm->dof, (long)spm->nnzexp, spm->layout );
967
968 switch(spm->fmttype){
969 case SpmCSC:
970 colsize = spm->n + 1;
971 rowsize = spm->nnz;
972 break;
973 case SpmCSR:
974 colsize = spm->nnz;
975 rowsize = spm->n + 1;
976 break;
977 case SpmIJV:
978 colsize = spm->nnz;
979 rowsize = spm->nnz;
980 break;
981 default:
982 colsize = 0;
983 rowsize = 0;
984 }
985
986 /*
987 * Write colptr
988 */
989 for (i=0; i<colsize; i++)
990 {
991 fprintf(outfile, "%ld ", (long)spm->colptr[i]);
992 if (i%4 == 3) fprintf(outfile, "\n");
993 }
994 if ((i-1)%4 !=3) fprintf(outfile, "\n");
995
996 /*
997 * Write rowptr
998 */
999 for (i=0; i<rowsize; i++)
1000 {
1001 fprintf(outfile, "%ld ", (long)spm->rowptr[i]);
1002 if (i%4 == 3) fprintf(outfile, "\n");
1003 }
1004 if ((i-1)%4 !=3) fprintf(outfile, "\n");
1005
1006 /*
1007 * Write dofs
1008 */
1009 if ( spm->dof <= 0 ) {
1010 for (i=0; i<spm->n+1; i++)
1011 {
1012 fprintf(outfile, "%ld ", (long)spm->dofs[i]);
1013 if (i%4 == 3) fprintf(outfile, "\n");
1014 }
1015 if ((i-1)%4 !=3) fprintf(outfile, "\n");
1016 }
1017
1018 /*
1019 * Write values
1020 */
1021 switch( spm->flttype ) {
1022 case SpmPattern:
1023 break;
1024 case SpmFloat:
1025 writeArrayOfFloat( outfile, spm->nnzexp, spm->values );
1026 break;
1027 case SpmDouble:
1028 writeArrayOfDouble( outfile, spm->nnzexp, spm->values );
1029 break;
1030 case SpmComplex32:
1031 writeArrayOfComplex32( outfile, spm->nnzexp, spm->values );
1032 break;
1033 case SpmComplex64:
1034 writeArrayOfComplex64( outfile, spm->nnzexp, spm->values );
1035 break;
1036 }
1037
1038 fclose(outfile);
1039 return SPM_SUCCESS;
1040}
1041
1042/**
1043 *******************************************************************************
1044 *
1045 * @ingroup spm
1046 *
1047 * @brief Save the spm structure into a file (internal format).
1048 *
1049 *******************************************************************************
1050 *
1051 * @param[in] spm
1052 * The sparse matrix to write into the file.
1053 *
1054 * @param[in] filename
1055 * The path where to store the spm. If filename == NULL, the spm
1056 * is saved into the "./matrix.spm" file.
1057 *
1058 ********************************************************************************
1059 *
1060 * @retval SPM_SUCCESS if the save happened successfully.
1061 * @retval SPM_ERR_FILE if the input format is incorrect.
1062 *
1063 *******************************************************************************/
1064int
1066 const char *filename )
1067{
1068 int rc = 0;
1069 int clustnum;
1070
1071 clustnum = spm->clustnum;
1072
1073 if ( clustnum == 0 ) {
1074 spmatrix_t *spm_local_ptr;
1075 spmatrix_t spm_local;
1076
1077#if defined(SPM_WITH_MPI)
1078 /* Gather the spm on one node */
1079 if( !spm->replicated ) {
1080 spmGather( spm, 0, &spm_local );
1081 spm_local_ptr = &spm_local;
1082 }
1083 else
1084#endif
1085 {
1086 spm_local_ptr = (spmatrix_t *)spm;
1087 }
1088 rc = spm_save_local( spm_local_ptr, filename );
1089
1090#if defined(SPM_WITH_MPI)
1091 MPI_Bcast( &rc, 1, MPI_INT, 0, spm->comm );
1092#endif
1093
1094 if ( !spm->replicated )
1095 {
1096 spmExit(&spm_local);
1097 }
1098 }
1099#if defined(SPM_WITH_MPI)
1100 else {
1101 assert( spm->clustnbr > 1 );
1102
1103 /* Gather the spm on one node */
1104 if( !spm->replicated ) {
1105 spmGather( spm, 0, NULL );
1106 }
1107 MPI_Bcast( &rc, 1, MPI_INT, 0, spm->comm );
1108 }
1109#endif
1110
1111 return rc;
1112}
static size_t spm_size_of(spm_coeftype_t type)
Double datatype that is not converted through precision generator functions.
Definition datatypes.h:209
float _Complex spm_complex32_t
Definition datatypes.h:161
enum spm_layout_e spm_layout_t
Direction of the matrix storage.
enum spm_fmttype_e spm_fmttype_t
Sparse matrix format.
enum spm_coeftype_e spm_coeftype_t
Arithmetic types.
enum spm_mtxtype_e spm_mtxtype_t
Matrix symmetry type property.
@ SPM_SUCCESS
Definition const.h:82
@ SPM_ERR_FILE
Definition const.h:90
@ SpmComplex64
Definition const.h:66
@ SpmDouble
Definition const.h:64
@ SpmFloat
Definition const.h:63
@ SpmComplex32
Definition const.h:65
@ SpmPattern
Definition const.h:62
@ SpmCSC
Definition const.h:73
@ SpmCSR
Definition const.h:74
@ SpmIJV
Definition const.h:75
spm_coeftype_t flttype
Definition spm.h:67
spm_int_t * dofs
Definition spm.h:85
void * values
Definition spm.h:92
int clustnbr
Definition spm.h:96
spm_int_t nnz
Definition spm.h:75
spm_int_t nnzexp
Definition spm.h:80
spm_int_t * rowptr
Definition spm.h:90
spm_fmttype_t fmttype
Definition spm.h:69
spm_int_t dof
Definition spm.h:82
spm_int_t n
Definition spm.h:73
spm_mtxtype_t mtxtype
Definition spm.h:65
spm_layout_t layout
Definition spm.h:87
SPM_Comm comm
Definition spm.h:97
spm_int_t * colptr
Definition spm.h:89
int clustnum
Definition spm.h:95
spm_int_t gN
Definition spm.h:72
int replicated
Definition spm.h:98
int spmGather(const spmatrix_t *spm_scattered, int root, spmatrix_t *opt_spm_gathered)
Gather a distributed Sparse Matrix on a single node.
Definition spm_gather.c:436
static int spm_save_local(const spmatrix_t *spm, const char *filename)
Save the global spm structure into a file (internal format).
Definition spm_io.c:942
int spmLoadDist(spmatrix_t *spm, const char *filename, SPM_Comm comm)
Load the spm structure from a file (internal format).
Definition spm_io.c:664
void spmExit(spmatrix_t *spm)
Cleanup the spm structure but do not free the spm pointer.
Definition spm.c:224
int spmLoad(spmatrix_t *spm, const char *filename)
Load the spm structure from a file (internal format).
Definition spm_io.c:749
void spmUpdateComputedFields(spmatrix_t *spm)
Update all the computed fields based on the static values stored.
int spmScatter(spmatrix_t *spm_scattered, int root, const spmatrix_t *opt_spm_gathered, spm_int_t opt_n, const spm_int_t *opt_loc2glob, int opt_distByColumn, SPM_Comm opt_comm)
Scatter the SPM thanks to loc2glob.
int spm_int_t
The main integer datatype used in spm arrays.
Definition datatypes.h:70
void spmInitDist(spmatrix_t *spm, SPM_Comm comm)
Init the spm structure with a specific communicator.
Definition spm.c:101
static int spm_load_local(spmatrix_t *spm, const char *filename)
Load the spm structure from a file (internal format). One node only.
Definition spm_io.c:464
int spmSave(const spmatrix_t *spm, const char *filename)
Save the spm structure into a file (internal format).
Definition spm_io.c:1065
The sparse matrix data structure.
Definition spm.h:64
static int readArrayOfComplex64(FILE *stream, spm_int_t n, spm_complex64_t *array)
Read an array of 64bits complex.
Definition spm_io.c:127
static int writeArrayOfFloat(FILE *outfile, spm_int_t n, const float *array)
write an array of float.
Definition spm_io.c:902
static int writeArrayOfDouble(FILE *outfile, spm_int_t n, const double *array)
write an array of double.
Definition spm_io.c:861
static int readArrayOfInteger(FILE *stream, spm_int_t n, spm_int_t *array)
Read an array of integer.
Definition spm_io.c:46
static int readArrayOfComplex32(FILE *stream, spm_int_t n, spm_complex32_t *array)
Read an array of 32bits complex.
Definition spm_io.c:213
static int writeArrayOfComplex64(FILE *outfile, spm_int_t n, const spm_complex64_t *array)
write an array of 64bits complex.
Definition spm_io.c:779
static int writeArrayOfComplex32(FILE *outfile, spm_int_t n, const spm_complex32_t *array)
write an array of 32bits complex.
Definition spm_io.c:820
static int readArrayOfDouble(FILE *stream, spm_int_t n, double *array)
Read an array of double.
Definition spm_io.c:299
static int readArrayOfFloat(FILE *stream, spm_int_t n, float *array)
Read an array of float.
Definition spm_io.c:383