50 long tmp1, tmp2, tmp3, tmp4;
54 for (i=0; i<(n-3); i+=4)
56 if (4 != fscanf(stream,
"%ld %ld %ld %ld", &tmp1, &tmp2, &tmp3, &tmp4)){
57 spm_print_error(
"spmLoad: Wrong input format");
71 if (3 != fscanf(stream,
"%ld %ld %ld", &tmp1, &tmp2, &tmp3)){
72 spm_print_error(
"spmLoad: Wrong input format");
81 if (2 != fscanf(stream,
"%ld %ld", &tmp1, &tmp2)){
82 spm_print_error(
"spmLoad: Wrong input format");
90 if (1 != fscanf(stream,
"%ld", &tmp1)){
91 spm_print_error(
"spmLoad: Wrong input format");
129 spm_complex64_t *array )
131 double tmp1, tmp2, tmp3, tmp4;
132 double tmp5, tmp6, tmp7, tmp8;
136 for (i=0; i<(n-3); i+=4)
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");
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);
154 if (6 != fscanf(stream,
"%lg %lg %lg %lg %lg %lg",
155 &tmp1, &tmp2, &tmp3, &tmp4,
157 spm_print_error(
"spmLoad: Wrong input format");
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);
166 if (4 != fscanf(stream,
"%lg %lg %lg %lg",
167 &tmp1, &tmp2, &tmp3, &tmp4)){
168 spm_print_error(
"spmLoad: Wrong input format");
171 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
172 array[i+1] = (spm_complex64_t)(tmp3 + I * tmp4);
176 if (2 != fscanf(stream,
"%lg %lg",
178 spm_print_error(
"spmLoad: Wrong input format");
181 array[i ] = (spm_complex64_t)(tmp1 + I * tmp2);
217 float tmp1, tmp2, tmp3, tmp4;
218 float tmp5, tmp6, tmp7, tmp8;
222 for (i=0; i<(n-3); i+=4)
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");
240 if (6 != fscanf(stream,
"%g %g %g %g %g %g",
241 &tmp1, &tmp2, &tmp3, &tmp4,
243 spm_print_error(
"spmLoad: Wrong input format");
252 if (4 != fscanf(stream,
"%g %g %g %g",
253 &tmp1, &tmp2, &tmp3, &tmp4)){
254 spm_print_error(
"spmLoad: Wrong input format");
262 if (2 != fscanf(stream,
"%g %g",
264 spm_print_error(
"spmLoad: Wrong input format");
303 double tmp1, tmp2, tmp3, tmp4;
307 for (i=0; i<(n-3); i+=4)
309 if (4 != fscanf(stream,
"%lg %lg %lg %lg",
310 &tmp1, &tmp2, &tmp3, &tmp4)){
311 spm_print_error(
"spmLoad: Wrong input format");
314 array[i ] = (double)(tmp1);
315 array[i+1] = (double)(tmp2);
316 array[i+2] = (double)(tmp3);
317 array[i+3] = (double)(tmp4);
324 if (1 != fscanf(stream,
"%lg %lg %lg",
325 &tmp1, &tmp2, &tmp3)){
326 spm_print_error(
"spmLoad: Wrong input format");
329 array[i ] = (double)(tmp1);
330 array[i+1] = (double)(tmp2);
331 array[i+2] = (double)(tmp3);
335 if (2 != fscanf(stream,
"%lg %lg",
337 spm_print_error(
"spmLoad: Wrong input format");
340 array[i ] = (double)(tmp1);
341 array[i+1] = (double)(tmp2);
345 if (1 != fscanf(stream,
"%lg",
347 spm_print_error(
"spmLoad: Wrong input format");
350 array[i ] = (double)(tmp1);
387 float tmp1, tmp2, tmp3, tmp4;
391 for (i=0; i<(n-3); i+=4)
393 if (4 != fscanf(stream,
"%g %g %g %g",
394 &tmp1, &tmp2, &tmp3, &tmp4)){
395 spm_print_error(
"spmLoad: Wrong input format");
398 array[i ] = (float)(tmp1);
399 array[i+1] = (float)(tmp2);
400 array[i+2] = (float)(tmp3);
401 array[i+3] = (float)(tmp4);
408 if (3 != fscanf(stream,
"%g %g %g",
409 &tmp1, &tmp2, &tmp3)){
410 spm_print_error(
"spmLoad: Wrong input format");
413 array[i ] = (float)(tmp1);
414 array[i+1] = (float)(tmp2);
415 array[i+2] = (float)(tmp3);
419 if (2 != fscanf(stream,
"%g %g",
421 spm_print_error(
"spmLoad: Wrong input format");
424 array[i ] = (float)(tmp1);
425 array[i+1] = (float)(tmp2);
429 if (1 != fscanf(stream,
"%g", &tmp1)){
430 spm_print_error(
"spmLoad: Wrong input format");
433 array[i ] = (float)(tmp1);
465 const char *filename )
475 if ( filename == NULL ) {
476 filename =
"matrix.spm";
479 infile = fopen( filename,
"r" );
480 if ( infile == NULL ) {
481 spm_print_error(
"spmLoad: Impossible to open the file matrix.spm\n");
489 test = fgets( line, 256, infile );
490 if ( test != line ) {
495 while( line[0] ==
'#' );
500 for (i=0; i<256; i++) {
501 if ( (line[i] == EOF ) || (line[i] ==
'\n') )
512 int version, mtxtype, flttype, fmttype, dof, layout;
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 ) )
523 if ( ( gN <= 0 ) || ( n <= 0 ) || ( nnz <= 0 ) )
548 colsize = spm->
n + 1;
553 rowsize = spm->
n + 1;
586 if ( spm->
dof > 0 ) {
603 assert( nnzexp == spm->
nnzexp );
604 assert( spm->
gN == gN );
665 const char *filename,
671 spm_only_with_mpi
int clustnbr = 1;
672 spm_only_with_mpi
int distbycol = -1;
674#if defined(SPM_WITH_MPI)
676 MPI_Comm_rank( comm, &clustnum );
677 MPI_Comm_size( comm, &clustnbr );
681 if ( clustnum == 0 ) {
687#if defined(SPM_WITH_MPI)
688 MPI_Bcast( &rc, 1, MPI_INT, 0, comm );
695#if defined(SPM_WITH_MPI)
703 (clustnum == 0) ? &spmlocal : NULL,
709 if ( clustnum == 0 ) {
717 assert( clustnum == 0 );
750 const char *filename )
752 return spmLoadDist( spm, filename, MPI_COMM_WORLD );
781 const spm_complex64_t *array )
788 fprintf(outfile,
"%e %e ", creal(array[i]), cimag(array[i]));
789 if (i%4 == 3) fprintf(outfile,
"\n");
791 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
829 fprintf(outfile,
"%e %e ", crealf(array[i]), cimagf(array[i]));
830 if (i%4 == 3) fprintf(outfile,
"\n");
832 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
863 const double *array )
870 fprintf(outfile,
"%e ", array[i]);
871 if (i%4 == 3) fprintf(outfile,
"\n");
873 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
911 fprintf(outfile,
"%e ", array[i]);
912 if (i%4 == 3) fprintf(outfile,
"\n");
914 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
943 const char *filename )
945 FILE *outfile = NULL;
948 if ( filename == NULL ) {
949 filename =
"matrix.spm";
952 outfile = fopen( filename,
"w" );
953 if ( outfile == NULL ) {
954 spm_print_error(
"spmSave: Impossible to open the file matrix.spm\n");
962 "# version mtxtype flttype fmttype gN n nnz dof nnzexp layout\n"
963 "%d %d %d %d %ld %ld %ld %d %ld %d\n",
965 (
long)spm->
gN, (
long)spm->
n, (
long)spm->
nnz,
970 colsize = spm->
n + 1;
975 rowsize = spm->
n + 1;
989 for (i=0; i<colsize; i++)
991 fprintf(outfile,
"%ld ", (
long)spm->
colptr[i]);
992 if (i%4 == 3) fprintf(outfile,
"\n");
994 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
999 for (i=0; i<rowsize; i++)
1001 fprintf(outfile,
"%ld ", (
long)spm->
rowptr[i]);
1002 if (i%4 == 3) fprintf(outfile,
"\n");
1004 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
1009 if ( spm->
dof <= 0 ) {
1010 for (i=0; i<spm->
n+1; i++)
1012 fprintf(outfile,
"%ld ", (
long)spm->
dofs[i]);
1013 if (i%4 == 3) fprintf(outfile,
"\n");
1015 if ((i-1)%4 !=3) fprintf(outfile,
"\n");
1066 const char *filename )
1073 if ( clustnum == 0 ) {
1077#if defined(SPM_WITH_MPI)
1081 spm_local_ptr = &spm_local;
1090#if defined(SPM_WITH_MPI)
1091 MPI_Bcast( &rc, 1, MPI_INT, 0, spm->
comm );
1099#if defined(SPM_WITH_MPI)
1107 MPI_Bcast( &rc, 1, MPI_INT, 0, spm->
comm );
static size_t spm_size_of(spm_coeftype_t type)
Double datatype that is not converted through precision generator functions.
float _Complex spm_complex32_t
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.
int spmGather(const spmatrix_t *spm_scattered, int root, spmatrix_t *opt_spm_gathered)
Gather a distributed Sparse Matrix on a single node.
static int spm_save_local(const spmatrix_t *spm, const char *filename)
Save the global spm structure into a file (internal format).
int spmLoadDist(spmatrix_t *spm, const char *filename, SPM_Comm comm)
Load the spm structure from a file (internal format).
void spmExit(spmatrix_t *spm)
Cleanup the spm structure but do not free the spm pointer.
int spmLoad(spmatrix_t *spm, const char *filename)
Load the spm structure from a file (internal format).
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.
void spmInitDist(spmatrix_t *spm, SPM_Comm comm)
Init the spm structure with a specific communicator.
static int spm_load_local(spmatrix_t *spm, const char *filename)
Load the spm structure from a file (internal format). One node only.
int spmSave(const spmatrix_t *spm, const char *filename)
Save the spm structure into a file (internal format).
The sparse matrix data structure.
static int readArrayOfComplex64(FILE *stream, spm_int_t n, spm_complex64_t *array)
Read an array of 64bits complex.
static int writeArrayOfFloat(FILE *outfile, spm_int_t n, const float *array)
write an array of float.
static int writeArrayOfDouble(FILE *outfile, spm_int_t n, const double *array)
write an array of double.
static int readArrayOfInteger(FILE *stream, spm_int_t n, spm_int_t *array)
Read an array of integer.
static int readArrayOfComplex32(FILE *stream, spm_int_t n, spm_complex32_t *array)
Read an array of 32bits complex.
static int writeArrayOfComplex64(FILE *outfile, spm_int_t n, const spm_complex64_t *array)
write an array of 64bits complex.
static int writeArrayOfComplex32(FILE *outfile, spm_int_t n, const spm_complex32_t *array)
write an array of 32bits complex.
static int readArrayOfDouble(FILE *stream, spm_int_t n, double *array)
Read an array of double.
static int readArrayOfFloat(FILE *stream, spm_int_t n, float *array)
Read an array of float.