-
Notifications
You must be signed in to change notification settings - Fork 0
/
hbio.c
1612 lines (1394 loc) · 59.9 KB
/
hbio.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* Harwell-Boeing File I/O in C
V. 1.0
National Institute of Standards and Technology, MD.
K.A. Remington
++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
NOTICE
Permission to use, copy, modify, and distribute this software and
its documentation for any purpose and without fee is hereby granted
provided that the above copyright notice appear in all copies and
that both the copyright notice and this permission notice appear in
supporting documentation.
Neither the Author nor the Institution (National Institute of Standards
and Technology) make any representations about the suitability of this
software for any purpose. This software is provided "as is" without
expressed or implied warranty.
++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
---------------------
INTERFACE DESCRIPTION
---------------------
---------------
QUERY FUNCTIONS
---------------
FUNCTION:
int readHB_info(const char *filename, int *M, int *N, int *nz,
char **Type, int *Nrhs)
DESCRIPTION:
The readHB_info function opens and reads the header information from
the specified Harwell-Boeing file, and reports back the number of rows
and columns in the stored matrix (M and N), the number of nonzeros in
the matrix (nz), the 3-character matrix type(Type), and the number of
right-hand-sides stored along with the matrix (Nrhs). This function
is designed to retrieve basic size information which can be used to
allocate arrays.
FUNCTION:
int readHB_header(FILE* in_file, char* Title, char* Key, char* Type,
int* Nrow, int* Ncol, int* Nnzero, int* Nrhs, int* Nrhsix,
char* Ptrfmt, char* Indfmt, char* Valfmt, char* Rhsfmt,
int* Ptrcrd, int* Indcrd, int* Valcrd, int* Rhscrd,
char *Rhstype)
DESCRIPTION:
More detailed than the readHB_info function, readHB_header() reads from
the specified Harwell-Boeing file all of the header information.
------------------------------
DOUBLE PRECISION I/O FUNCTIONS
------------------------------
FUNCTION:
int readHB_newmat_double(const char *filename, int *M, int *N, *int nz,
int **colptr, int **rowind, double**val)
int readHB_mat_double(const char *filename, int *colptr, int *rowind,
double*val)
DESCRIPTION:
This function opens and reads the specified file, interpreting its
contents as a sparse matrix stored in the Harwell/Boeing standard
format. (See readHB_aux_double to read auxillary vectors.)
-- Values are interpreted as double precision numbers. --
The "mat" function uses _pre-allocated_ vectors to hold the index and
nonzero value information.
The "newmat" function allocates vectors to hold the index and nonzero
value information, and returns pointers to these vectors along with
matrix dimension and number of nonzeros.
FUNCTION:
int readHB_aux_double(const char* filename, const char AuxType, double b[])
int readHB_newaux_double(const char* filename, const char AuxType, double** b)
DESCRIPTION:
This function opens and reads from the specified file auxillary vector(s).
The char argument Auxtype determines which type of auxillary vector(s)
will be read (if present in the file).
AuxType = 'F' right-hand-side
AuxType = 'G' initial estimate (Guess)
AuxType = 'X' eXact solution
If Nrhs > 1, all of the Nrhs vectors of the given type are read and
stored in column-major order in the vector b.
The "newaux" function allocates a vector to hold the values retrieved.
The "mat" function uses a _pre-allocated_ vector to hold the values.
FUNCTION:
int writeHB_mat_double(const char* filename, int M, int N,
int nz, const int colptr[], const int rowind[],
const double val[], int Nrhs, const double rhs[],
const double guess[], const double exact[],
const char* Title, const char* Key, const char* Type,
char* Ptrfmt, char* Indfmt, char* Valfmt, char* Rhsfmt,
const char* Rhstype)
DESCRIPTION:
The writeHB_mat_double function opens the named file and writes the specified
matrix and optional auxillary vector(s) to that file in Harwell-Boeing
format. The format arguments (Ptrfmt,Indfmt,Valfmt, and Rhsfmt) are
character strings specifying "Fortran-style" output formats -- as they
would appear in a Harwell-Boeing file. They are used to produce output
which is as close as possible to what would be produced by Fortran code,
but note that "D" and "P" edit descriptors are not supported.
If NULL, the following defaults will be used:
Ptrfmt = Indfmt = "(8I10)"
Valfmt = Rhsfmt = "(4E20.13)"
-----------------------
CHARACTER I/O FUNCTIONS
-----------------------
FUNCTION:
int readHB_mat_char(const char* filename, int colptr[], int rowind[],
char val[], char* Valfmt)
int readHB_newmat_char(const char* filename, int* M, int* N, int* nonzeros,
int** colptr, int** rowind, char** val, char** Valfmt)
DESCRIPTION:
This function opens and reads the specified file, interpreting its
contents as a sparse matrix stored in the Harwell/Boeing standard
format. (See readHB_aux_char to read auxillary vectors.)
-- Values are interpreted as char strings. --
(Used to translate exact values from the file into a new storage format.)
The "mat" function uses _pre-allocated_ arrays to hold the index and
nonzero value information.
The "newmat" function allocates char arrays to hold the index
and nonzero value information, and returns pointers to these arrays
along with matrix dimension and number of nonzeros.
FUNCTION:
int readHB_aux_char(const char* filename, const char AuxType, char b[])
int readHB_newaux_char(const char* filename, const char AuxType, char** b,
char** Rhsfmt)
DESCRIPTION:
This function opens and reads from the specified file auxillary vector(s).
The char argument Auxtype determines which type of auxillary vector(s)
will be read (if present in the file).
AuxType = 'F' right-hand-side
AuxType = 'G' initial estimate (Guess)
AuxType = 'X' eXact solution
If Nrhs > 1, all of the Nrhs vectors of the given type are read and
stored in column-major order in the vector b.
The "newaux" function allocates a character array to hold the values
retrieved.
The "mat" function uses a _pre-allocated_ array to hold the values.
FUNCTION:
int writeHB_mat_char(const char* filename, int M, int N,
int nz, const int colptr[], const int rowind[],
const char val[], int Nrhs, const char rhs[],
const char guess[], const char exact[],
const char* Title, const char* Key, const char* Type,
char* Ptrfmt, char* Indfmt, char* Valfmt, char* Rhsfmt,
const char* Rhstype)
DESCRIPTION:
The writeHB_mat_char function opens the named file and writes the specified
matrix and optional auxillary vector(s) to that file in Harwell-Boeing
format. The format arguments (Ptrfmt,Indfmt,Valfmt, and Rhsfmt) are
character strings specifying "Fortran-style" output formats -- as they
would appear in a Harwell-Boeing file. Valfmt and Rhsfmt must accurately
represent the character representation of the values stored in val[]
and rhs[].
If NULL, the following defaults will be used for the integer vectors:
Ptrfmt = Indfmt = "(8I10)"
Valfmt = Rhsfmt = "(4E20.13)"
*/
/*---------------------------------------------------------------------*/
/* If zero-based indexing is desired, _SP_base should be set to 0 */
/* This will cause indices read from H-B files to be decremented by 1 */
/* and indices written to H-B files to be incremented by 1 */
/* <<< Standard usage is _SP_base = 1 >>> */
#ifndef _SP_base
#define _SP_base 1
#endif
/*---------------------------------------------------------------------*/
#include "hbio.h"
#include<stdio.h>
#include<stdlib.h>
#include<string.h>
#include<math.h>
#include<sys/types.h>
#include<malloc.h>
char *substr(const char* S, const int pos, const int len);
void upcase(char* S);
void IOHBTerminate(char* message);
int readHB_info(const char* filename, int* M, int* N, int* nz, char** Type,
int* Nrhs)
{
/****************************************************************************/
/* The readHB_info function opens and reads the header information from */
/* the specified Harwell-Boeing file, and reports back the number of rows */
/* and columns in the stored matrix (M and N), the number of nonzeros in */
/* the matrix (nz), and the number of right-hand-sides stored along with */
/* the matrix (Nrhs). */
/* */
/* For a description of the Harwell Boeing standard, see: */
/* Duff, et al., ACM TOMS Vol.15, No.1, March 1989 */
/* */
/* ---------- */
/* **CAVEAT** */
/* ---------- */
/* ** If the input file does not adhere to the H/B format, the ** */
/* ** results will be unpredictable. ** */
/* */
/****************************************************************************/
FILE *in_file;
int Ptrcrd, Indcrd, Valcrd, Rhscrd;
int Nrow, Ncol, Nnzero, Nrhsix;
char *mat_type;
char Title[73], Key[9], Rhstype[4];
char Ptrfmt[17], Indfmt[17], Valfmt[21], Rhsfmt[21];
mat_type = (char *) malloc(4);
if ( mat_type == NULL ) IOHBTerminate("Insufficient memory for mat_typen");
if ( (in_file = fopen( filename, "r")) == NULL ) {
fprintf(stderr,"Error: Cannot open file: %s\n",filename);
return 0;
}
readHB_header(in_file, Title, Key, mat_type, &Nrow, &Ncol, &Nnzero,
Nrhs, &Nrhsix,
Ptrfmt, Indfmt, Valfmt, Rhsfmt,
&Ptrcrd, &Indcrd, &Valcrd, &Rhscrd, Rhstype);
fclose(in_file);
*Type = mat_type;
*(*Type+3) = (char) NULL;
*M = Nrow;
*N = Ncol;
*nz = Nnzero;
if (Rhscrd == 0) {*Nrhs = 0;}
/* In verbose mode, print some of the header information: */
/*
if (verbose == 1)
{
printf("Reading from Harwell-Boeing file %s (verbose on)...\n",filename);
printf(" Title: %s\n",Title);
printf(" Key: %s\n",Key);
printf(" The stored matrix is %i by %i with %i nonzeros.\n",
*M, *N, *nz );
printf(" %i right-hand--side(s) stored.\n",*Nrhs);
}
*/
return 1;
}
int readHB_header(FILE* in_file, char* Title, char* Key, char* Type,
int* Nrow, int* Ncol, int* Nnzero, int* Nrhs, int* Nrhsix,
char* Ptrfmt, char* Indfmt, char* Valfmt, char* Rhsfmt,
int* Ptrcrd, int* Indcrd, int* Valcrd, int* Rhscrd,
char *Rhstype)
{
/*************************************************************************/
/* Read header information from the named H/B file... */
/*************************************************************************/
int Totcrd,Neltvl;
char line[BUFSIZ];
/* First line: */
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) first line of HB file.\n");
(void) sscanf(line, "%72c%8[^\n]", Title, Key);
*(Key+8) = (char) NULL;
*(Title+72) = (char) NULL;
/* Second line: */
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) second line of HB file.\n");
if ( sscanf(line,"%i",&Totcrd) != 1) Totcrd = 0;
if ( sscanf(line,"%*i%i",Ptrcrd) != 1) *Ptrcrd = 0;
if ( sscanf(line,"%*i%*i%i",Indcrd) != 1) *Indcrd = 0;
if ( sscanf(line,"%*i%*i%*i%i",Valcrd) != 1) *Valcrd = 0;
if ( sscanf(line,"%*i%*i%*i%*i%i",Rhscrd) != 1) *Rhscrd = 0;
/* Third line: */
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) third line of HB file.\n");
if ( sscanf(line, "%3c", Type) != 1)
IOHBTerminate("iohb.c: Invalid Type info, line 3 of Harwell-Boeing file.\n");
upcase(Type);
if ( sscanf(line,"%*3c%i",Nrow) != 1) *Nrow = 0 ;
if ( sscanf(line,"%*3c%*i%i",Ncol) != 1) *Ncol = 0 ;
if ( sscanf(line,"%*3c%*i%*i%i",Nnzero) != 1) *Nnzero = 0 ;
if ( sscanf(line,"%*3c%*i%*i%*i%i",&Neltvl) != 1) Neltvl = 0 ;
/* Fourth line: */
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) fourth line of HB file.\n");
if ( sscanf(line, "%16c",Ptrfmt) != 1)
IOHBTerminate("iohb.c: Invalid format info, line 4 of Harwell-Boeing file.\n");
if ( sscanf(line, "%*16c%16c",Indfmt) != 1)
IOHBTerminate("iohb.c: Invalid format info, line 4 of Harwell-Boeing file.\n");
if ( sscanf(line, "%*16c%*16c%20c",Valfmt) != 1)
IOHBTerminate("iohb.c: Invalid format info, line 4 of Harwell-Boeing file.\n");
sscanf(line, "%*16c%*16c%*20c%20c",Rhsfmt);
*(Ptrfmt+16) = (char) NULL;
*(Indfmt+16) = (char) NULL;
*(Valfmt+20) = (char) NULL;
*(Rhsfmt+20) = (char) NULL;
/* (Optional) Fifth line: */
if (*Rhscrd != 0 )
{
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) fifth line of HB file.\n");
if ( sscanf(line, "%3c", Rhstype) != 1)
IOHBTerminate("iohb.c: Invalid RHS type information, line 5 of Harwell-Boeing file.\n");
if ( sscanf(line, "%*3c%i", Nrhs) != 1) *Nrhs = 0;
if ( sscanf(line, "%*3c%*i%i", Nrhsix) != 1) *Nrhsix = 0;
}
return 1;
}
int readHB_mat_double(const char* filename, int colptr[], int rowind[],
double val[])
{
/****************************************************************************/
/* This function opens and reads the specified file, interpreting its */
/* contents as a sparse matrix stored in the Harwell/Boeing standard */
/* format and creating compressed column storage scheme vectors to hold */
/* the index and nonzero value information. */
/* */
/* ---------- */
/* **CAVEAT** */
/* ---------- */
/* Parsing real formats from Fortran is tricky, and this file reader */
/* does not claim to be foolproof. It has been tested for cases when */
/* the real values are printed consistently and evenly spaced on each */
/* line, with Fixed (F), and Exponential (E or D) formats. */
/* */
/* ** If the input file does not adhere to the H/B format, the ** */
/* ** results will be unpredictable. ** */
/* */
/****************************************************************************/
FILE *in_file;
int i,j,ind,col,offset,count,last,Nrhs,Nrhsix;
int Ptrcrd, Indcrd, Valcrd, Rhscrd;
int Nrow, Ncol, Nnzero, Nentries;
int Ptrperline, Ptrwidth, Indperline, Indwidth;
int Valperline, Valwidth, Valprec;
int Valflag; /* Indicates 'E','D', or 'F' float format */
char* ThisElement;
char Title[73], Key[9], Type[4], Rhstype[4];
char Ptrfmt[17], Indfmt[17], Valfmt[21], Rhsfmt[21];
char line[BUFSIZ];
if ( (in_file = fopen( filename, "r")) == NULL ) {
fprintf(stderr,"Error: Cannot open file: %s\n",filename);
return 0;
}
readHB_header(in_file, Title, Key, Type, &Nrow, &Ncol, &Nnzero,
&Nrhs, &Nrhsix,
Ptrfmt, Indfmt, Valfmt, Rhsfmt,
&Ptrcrd, &Indcrd, &Valcrd, &Rhscrd, Rhstype);
/* Parse the array input formats from Line 3 of HB file */
ParseIfmt(Ptrfmt,&Ptrperline,&Ptrwidth);
ParseIfmt(Indfmt,&Indperline,&Indwidth);
if ( Type[0] != 'P' ) { /* Skip if pattern only */
ParseRfmt(Valfmt,&Valperline,&Valwidth,&Valprec,&Valflag);
}
/* Read column pointer array: */
offset = 1-_SP_base; /* if base 0 storage is declared (via macro definition), */
/* then storage entries are offset by 1 */
ThisElement = (char *) malloc(Ptrwidth+1);
if ( ThisElement == NULL ) IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Ptrwidth) = (char) NULL;
count=0;
for (i=0;i<Ptrcrd;i++)
{
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) line in pointer data region of HB file.\n");
col = 0;
for (ind = 0;ind<Ptrperline;ind++)
{
if (count > Ncol) break;
strncpy(ThisElement,line+col,Ptrwidth);
/* ThisElement = substr(line,col,Ptrwidth); */
colptr[count] = atoi(ThisElement)-offset;
count++; col += Ptrwidth;
}
}
free(ThisElement);
/* Read row index array: */
ThisElement = (char *) malloc(Indwidth+1);
if ( ThisElement == NULL ) IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Indwidth) = (char) NULL;
count = 0;
for (i=0;i<Indcrd;i++)
{
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) line in index data region of HB file.\n");
col = 0;
for (ind = 0;ind<Indperline;ind++)
{
if (count == Nnzero) break;
strncpy(ThisElement,line+col,Indwidth);
/* ThisElement = substr(line,col,Indwidth); */
rowind[count] = atoi(ThisElement)-offset;
count++; col += Indwidth;
}
}
free(ThisElement);
/* Read array of values: */
if ( Type[0] != 'P' ) { /* Skip if pattern only */
if ( Type[0] == 'C' ) Nentries = 2*Nnzero;
else Nentries = Nnzero;
ThisElement = (char *) malloc(Valwidth+1);
if ( ThisElement == NULL )
IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Valwidth) = (char) NULL;
count = 0;
for (i=0;i<Valcrd;i++)
{
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) line in value data region of HB file.\n");
if (Valflag == 'D') {
while( strchr(line,'D') ) *strchr(line,'D') = 'E';
/* *strchr(Valfmt,'D') = 'E'; */
}
col = 0;
for (ind = 0;ind<Valperline;ind++)
{
if (count == Nentries) break;
strncpy(ThisElement,line+col,Valwidth);
/*ThisElement = substr(line,col,Valwidth);*/
if ( Valflag != 'F' && strchr(ThisElement,'E') == NULL ) {
/* insert a char prefix for exp */
last = strlen(ThisElement);
for (j=last+1;j>=0;j--) {
ThisElement[j] = ThisElement[j-1];
if ( ThisElement[j] == '+' || ThisElement[j] == '-' ) {
ThisElement[j-1] = Valflag;
break;
}
}
}
val[count] = atof(ThisElement);
count++; col += Valwidth;
}
}
free(ThisElement);
}
fclose(in_file);
return 1;
}
int readHB_newmat_double(const char* filename, int* M, int* N, int* nonzeros,
int** colptr, int** rowind, double** val)
{
int Nrhs;
char *Type;
readHB_info(filename, M, N, nonzeros, &Type, &Nrhs);
*colptr = (int *)malloc((*N+1)*sizeof(int));
if ( *colptr == NULL ) IOHBTerminate("Insufficient memory for colptr.\n");
*rowind = (int *)malloc(*nonzeros*sizeof(int));
if ( *rowind == NULL ) IOHBTerminate("Insufficient memory for rowind.\n");
if ( Type[0] == 'C' ) {
/*
fprintf(stderr, "Warning: Reading complex data from HB file %s.\n",filename);
fprintf(stderr, " Real and imaginary parts will be interlaced in val[].\n");
*/
/* Malloc enough space for real AND imaginary parts of val[] */
*val = (double *)malloc(*nonzeros*sizeof(double)*2);
if ( *val == NULL ) IOHBTerminate("Insufficient memory for val.\n");
} else {
if ( Type[0] != 'P' ) {
/* Malloc enough space for real array val[] */
*val = (double *)malloc(*nonzeros*sizeof(double));
if ( *val == NULL ) IOHBTerminate("Insufficient memory for val.\n");
}
} /* No val[] space needed if pattern only */
return readHB_mat_double(filename, *colptr, *rowind, *val);
}
int readHB_aux_double(const char* filename, const char AuxType, double b[])
{
/****************************************************************************/
/* This function opens and reads the specified file, placing auxillary */
/* vector(s) of the given type (if available) in b. */
/* Return value is the number of vectors successfully read. */
/* */
/* AuxType = 'F' full right-hand-side vector(s) */
/* AuxType = 'G' initial Guess vector(s) */
/* AuxType = 'X' eXact solution vector(s) */
/* */
/* ---------- */
/* **CAVEAT** */
/* ---------- */
/* Parsing real formats from Fortran is tricky, and this file reader */
/* does not claim to be foolproof. It has been tested for cases when */
/* the real values are printed consistently and evenly spaced on each */
/* line, with Fixed (F), and Exponential (E or D) formats. */
/* */
/* ** If the input file does not adhere to the H/B format, the ** */
/* ** results will be unpredictable. ** */
/* */
/****************************************************************************/
FILE *in_file;
int i,j,n,maxcol,start,stride,col,last,linel;
int Ptrcrd, Indcrd, Valcrd, Rhscrd;
int Nrow, Ncol, Nnzero, Nentries;
int Nrhs, Nrhsix, nvecs, rhsi;
int Rhsperline, Rhswidth, Rhsprec;
int Rhsflag;
char *ThisElement;
char Title[73], Key[9], Type[4], Rhstype[4];
char Ptrfmt[17], Indfmt[17], Valfmt[21], Rhsfmt[21];
char line[BUFSIZ];
if ((in_file = fopen( filename, "r")) == NULL) {
fprintf(stderr,"Error: Cannot open file: %s\n",filename);
return 0;
}
readHB_header(in_file, Title, Key, Type, &Nrow, &Ncol, &Nnzero,
&Nrhs, &Nrhsix,
Ptrfmt, Indfmt, Valfmt, Rhsfmt,
&Ptrcrd, &Indcrd, &Valcrd, &Rhscrd, Rhstype);
if (Nrhs <= 0)
{
fprintf(stderr, "Warn: Attempt to read auxillary vector(s) when none are present.\n");
return 0;
}
if (Rhstype[0] != 'F' )
{
fprintf(stderr,"Warn: Attempt to read auxillary vector(s) which are not stored in Full form.\n");
fprintf(stderr," Rhs must be specified as full. \n");
return 0;
}
/* If reading complex data, allow for interleaved real and imaginary values. */
if ( Type[0] == 'C' ) {
Nentries = 2*Nrow;
} else {
Nentries = Nrow;
}
nvecs = 1;
if ( Rhstype[1] == 'G' ) nvecs++;
if ( Rhstype[2] == 'X' ) nvecs++;
if ( AuxType == 'G' && Rhstype[1] != 'G' ) {
fprintf(stderr, "Warn: Attempt to read auxillary Guess vector(s) when none are present.\n");
return 0;
}
if ( AuxType == 'X' && Rhstype[2] != 'X' ) {
fprintf(stderr, "Warn: Attempt to read auxillary eXact solution vector(s) when none are present.\n");
return 0;
}
ParseRfmt(Rhsfmt, &Rhsperline, &Rhswidth, &Rhsprec,&Rhsflag);
maxcol = Rhsperline*Rhswidth;
/* Lines to skip before starting to read RHS values... */
n = Ptrcrd + Indcrd + Valcrd;
for (i = 0; i < n; i++)
fgets(line, BUFSIZ, in_file);
/* start - number of initial aux vector entries to skip */
/* to reach first vector requested */
/* stride - number of aux vector entries to skip between */
/* requested vectors */
if ( AuxType == 'F' ) start = 0;
else if ( AuxType == 'G' ) start = Nentries;
else start = (nvecs-1)*Nentries;
stride = (nvecs-1)*Nentries;
fgets(line, BUFSIZ, in_file);
linel= strchr(line,'\n')-line;
col = 0;
/* Skip to initial offset */
for (i=0;i<start;i++) {
if ( col >= ( maxcol<linel?maxcol:linel ) ) {
fgets(line, BUFSIZ, in_file);
linel= strchr(line,'\n')-line;
col = 0;
}
col += Rhswidth;
}
if (Rhsflag == 'D') {
while( strchr(line,'D') ) *strchr(line,'D') = 'E';
}
/* Read a vector of desired type, then skip to next */
/* repeating to fill Nrhs vectors */
ThisElement = (char *) malloc(Rhswidth+1);
if ( ThisElement == NULL ) IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Rhswidth) = (char) NULL;
for (rhsi=0;rhsi<Nrhs;rhsi++) {
for (i=0;i<Nentries;i++) {
if ( col >= ( maxcol<linel?maxcol:linel ) ) {
fgets(line, BUFSIZ, in_file);
linel= strchr(line,'\n')-line;
if (Rhsflag == 'D') {
while( strchr(line,'D') ) *strchr(line,'D') = 'E';
}
col = 0;
}
strncpy(ThisElement,line+col,Rhswidth);
/*ThisElement = substr(line, col, Rhswidth);*/
if ( Rhsflag != 'F' && strchr(ThisElement,'E') == NULL ) {
/* insert a char prefix for exp */
last = strlen(ThisElement);
for (j=last+1;j>=0;j--) {
ThisElement[j] = ThisElement[j-1];
if ( ThisElement[j] == '+' || ThisElement[j] == '-' ) {
ThisElement[j-1] = Rhsflag;
break;
}
}
}
b[i] = atof(ThisElement);
col += Rhswidth;
}
/* Skip any interleaved Guess/eXact vectors */
for (i=0;i<stride;i++) {
if ( col >= ( maxcol<linel?maxcol:linel ) ) {
fgets(line, BUFSIZ, in_file);
linel= strchr(line,'\n')-line;
col = 0;
}
col += Rhswidth;
}
}
free(ThisElement);
fclose(in_file);
return Nrhs;
}
int readHB_newaux_double(const char* filename, const char AuxType, double** b)
{
int Nrhs,M,N,nonzeros;
char *Type;
readHB_info(filename, &M, &N, &nonzeros, &Type, &Nrhs);
if ( Nrhs <= 0 ) {
fprintf(stderr,"Warn: Requested read of aux vector(s) when none are present.\n");
return 0;
} else {
if ( Type[0] == 'C' ) {
fprintf(stderr, "Warning: Reading complex aux vector(s) from HB file %s.",filename);
fprintf(stderr, " Real and imaginary parts will be interlaced in b[].");
*b = (double *)malloc(M*Nrhs*sizeof(double)*2);
if ( *b == NULL ) IOHBTerminate("Insufficient memory for rhs.\n");
return readHB_aux_double(filename, AuxType, *b);
} else {
*b = (double *)malloc(M*Nrhs*sizeof(double));
if ( *b == NULL ) IOHBTerminate("Insufficient memory for rhs.\n");
return readHB_aux_double(filename, AuxType, *b);
}
}
}
int writeHB_mat_double(const char* filename, int M, int N,
int nz, const int colptr[], const int rowind[],
const double val[], int Nrhs, const double rhs[],
const double guess[], const double exact[],
const char* Title, const char* Key, const char* Type,
char* Ptrfmt, char* Indfmt, char* Valfmt, char* Rhsfmt,
const char* Rhstype)
{
/****************************************************************************/
/* The writeHB function opens the named file and writes the specified */
/* matrix and optional right-hand-side(s) to that file in Harwell-Boeing */
/* format. */
/* */
/* For a description of the Harwell Boeing standard, see: */
/* Duff, et al., ACM TOMS Vol.15, No.1, March 1989 */
/* */
/****************************************************************************/
FILE *out_file;
int i,j,entry,offset,acount,linemod;
int totcrd, ptrcrd, indcrd, valcrd, rhscrd;
int nvalentries, nrhsentries;
int Ptrperline, Ptrwidth, Indperline, Indwidth;
int Rhsperline, Rhswidth, Rhsprec;
int Rhsflag;
int Valperline, Valwidth, Valprec;
int Valflag; /* Indicates 'E','D', or 'F' float format */
char pformat[16],iformat[16],vformat[19],rformat[19];
if ( Type[0] == 'C' ) {
nvalentries = 2*nz;
nrhsentries = 2*M;
} else {
nvalentries = nz;
nrhsentries = M;
}
if ( filename != NULL ) {
if ( (out_file = fopen( filename, "w")) == NULL ) {
fprintf(stderr,"Error: Cannot open file: %s\n",filename);
return 0;
}
} else out_file = stdout;
if ( Ptrfmt == NULL ) Ptrfmt = "(8I10)";
ParseIfmt(Ptrfmt,&Ptrperline,&Ptrwidth);
sprintf(pformat,"%%%dd",Ptrwidth);
ptrcrd = (N+1)/Ptrperline;
if ( (N+1)%Ptrperline != 0) ptrcrd++;
if ( Indfmt == NULL ) Indfmt = Ptrfmt;
ParseIfmt(Indfmt,&Indperline,&Indwidth);
sprintf(iformat,"%%%dd",Indwidth);
indcrd = nz/Indperline;
if ( nz%Indperline != 0) indcrd++;
if ( Type[0] != 'P' ) { /* Skip if pattern only */
if ( Valfmt == NULL ) Valfmt = "(4E20.13)";
ParseRfmt(Valfmt,&Valperline,&Valwidth,&Valprec,&Valflag);
if (Valflag == 'D') *strchr(Valfmt,'D') = 'E';
if (Valflag == 'F')
sprintf(vformat,"%% %d.%df",Valwidth,Valprec);
else
sprintf(vformat,"%% %d.%dE",Valwidth,Valprec);
valcrd = nvalentries/Valperline;
if ( nvalentries%Valperline != 0) valcrd++;
} else valcrd = 0;
if ( Nrhs > 0 ) {
if ( Rhsfmt == NULL ) Rhsfmt = Valfmt;
ParseRfmt(Rhsfmt,&Rhsperline,&Rhswidth,&Rhsprec, &Rhsflag);
if (Rhsflag == 'F')
sprintf(rformat,"%% %d.%df",Rhswidth,Rhsprec);
else
sprintf(rformat,"%% %d.%dE",Rhswidth,Rhsprec);
if (Rhsflag == 'D') *strchr(Rhsfmt,'D') = 'E';
rhscrd = nrhsentries/Rhsperline;
if ( nrhsentries%Rhsperline != 0) rhscrd++;
if ( Rhstype[1] == 'G' ) rhscrd+=rhscrd;
if ( Rhstype[2] == 'X' ) rhscrd+=rhscrd;
rhscrd*=Nrhs;
} else rhscrd = 0;
totcrd = 4+ptrcrd+indcrd+valcrd+rhscrd;
/* Print header information: */
fprintf(out_file,"%-72s%-8s\n%14d%14d%14d%14d%14d\n",Title, Key, totcrd,
ptrcrd, indcrd, valcrd, rhscrd);
fprintf(out_file,"%3s%11s%14d%14d%14d\n",Type," ", M, N, nz);
fprintf(out_file,"%-16s%-16s%-20s", Ptrfmt, Indfmt, Valfmt);
if ( Nrhs != 0 ) {
/* Print Rhsfmt on fourth line and */
/* optional fifth header line for auxillary vector information: */
fprintf(out_file,"%-20s\n%-14s%d\n",Rhsfmt,Rhstype,Nrhs);
} else fprintf(out_file,"\n");
offset = 1-_SP_base; /* if base 0 storage is declared (via macro definition), */
/* then storage entries are offset by 1 */
/* Print column pointers: */
for (i=0;i<N+1;i++)
{
entry = colptr[i]+offset;
fprintf(out_file,pformat,entry);
if ( (i+1)%Ptrperline == 0 ) fprintf(out_file,"\n");
}
if ( (N+1) % Ptrperline != 0 ) fprintf(out_file,"\n");
/* Print row indices: */
for (i=0;i<nz;i++)
{
entry = rowind[i]+offset;
fprintf(out_file,iformat,entry);
if ( (i+1)%Indperline == 0 ) fprintf(out_file,"\n");
}
if ( nz % Indperline != 0 ) fprintf(out_file,"\n");
/* Print values: */
if ( Type[0] != 'P' ) { /* Skip if pattern only */
for (i=0;i<nvalentries;i++)
{
fprintf(out_file,vformat,val[i]);
if ( (i+1)%Valperline == 0 ) fprintf(out_file,"\n");
}
if ( nvalentries % Valperline != 0 ) fprintf(out_file,"\n");
/* If available, print right hand sides,
guess vectors and exact solution vectors: */
acount = 1;
linemod = 0;
if ( Nrhs > 0 ) {
for (i=0;i<Nrhs;i++)
{
for ( j=0;j<nrhsentries;j++ ) {
fprintf(out_file,rformat,rhs[j]);
if ( acount++%Rhsperline == linemod ) fprintf(out_file,"\n");
}
if ( acount%Rhsperline != linemod ) {
fprintf(out_file,"\n");
linemod = (acount-1)%Rhsperline;
}
rhs += nrhsentries;
if ( Rhstype[1] == 'G' ) {
for ( j=0;j<nrhsentries;j++ ) {
fprintf(out_file,rformat,guess[j]);
if ( acount++%Rhsperline == linemod ) fprintf(out_file,"\n");
}
if ( acount%Rhsperline != linemod ) {
fprintf(out_file,"\n");
linemod = (acount-1)%Rhsperline;
}
guess += nrhsentries;
}
if ( Rhstype[2] == 'X' ) {
for ( j=0;j<nrhsentries;j++ ) {
fprintf(out_file,rformat,exact[j]);
if ( acount++%Rhsperline == linemod ) fprintf(out_file,"\n");
}
if ( acount%Rhsperline != linemod ) {
fprintf(out_file,"\n");
linemod = (acount-1)%Rhsperline;
}
exact += nrhsentries;
}
}
}
}
if ( fclose(out_file) != 0){
fprintf(stderr,"Error closing file in writeHB_mat_double().\n");
return 0;
} else return 1;
}
int readHB_mat_char(const char* filename, int colptr[], int rowind[],
char val[], char* Valfmt)
{
/****************************************************************************/
/* This function opens and reads the specified file, interpreting its */
/* contents as a sparse matrix stored in the Harwell/Boeing standard */
/* format and creating compressed column storage scheme vectors to hold */
/* the index and nonzero value information. */
/* */
/* ---------- */
/* **CAVEAT** */
/* ---------- */
/* Parsing real formats from Fortran is tricky, and this file reader */
/* does not claim to be foolproof. It has been tested for cases when */
/* the real values are printed consistently and evenly spaced on each */
/* line, with Fixed (F), and Exponential (E or D) formats. */
/* */
/* ** If the input file does not adhere to the H/B format, the ** */
/* ** results will be unpredictable. ** */
/* */
/****************************************************************************/
FILE *in_file;
int i,j,ind,col,offset,count,last;
int Nrow,Ncol,Nnzero,Nentries,Nrhs,Nrhsix;
int Ptrcrd, Indcrd, Valcrd, Rhscrd;
int Ptrperline, Ptrwidth, Indperline, Indwidth;
int Valperline, Valwidth, Valprec;
int Valflag; /* Indicates 'E','D', or 'F' float format */
char* ThisElement;
char line[BUFSIZ];
char Title[73], Key[9], Type[4], Rhstype[4];
char Ptrfmt[17], Indfmt[17], Rhsfmt[21];
if ( (in_file = fopen( filename, "r")) == NULL ) {
fprintf(stderr,"Error: Cannot open file: %s\n",filename);
return 0;
}
readHB_header(in_file, Title, Key, Type, &Nrow, &Ncol, &Nnzero,
&Nrhs, &Nrhsix,
Ptrfmt, Indfmt, Valfmt, Rhsfmt,
&Ptrcrd, &Indcrd, &Valcrd, &Rhscrd, Rhstype);
/* Parse the array input formats from Line 3 of HB file */
ParseIfmt(Ptrfmt,&Ptrperline,&Ptrwidth);
ParseIfmt(Indfmt,&Indperline,&Indwidth);
if ( Type[0] != 'P' ) { /* Skip if pattern only */
ParseRfmt(Valfmt,&Valperline,&Valwidth,&Valprec,&Valflag);
if (Valflag == 'D') {
*strchr(Valfmt,'D') = 'E';
}
}
/* Read column pointer array: */
offset = 1-_SP_base; /* if base 0 storage is declared (via macro definition), */
/* then storage entries are offset by 1 */
ThisElement = (char *) malloc(Ptrwidth+1);
if ( ThisElement == NULL ) IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Ptrwidth) = (char) NULL;
count=0;
for (i=0;i<Ptrcrd;i++)
{
fgets(line, BUFSIZ, in_file);
if ( sscanf(line,"%*s") < 0 )
IOHBTerminate("iohb.c: Null (or blank) line in pointer data region of HB file.\n");
col = 0;
for (ind = 0;ind<Ptrperline;ind++)
{
if (count > Ncol) break;
strncpy(ThisElement,line+col,Ptrwidth);
/*ThisElement = substr(line,col,Ptrwidth);*/
colptr[count] = atoi(ThisElement)-offset;
count++; col += Ptrwidth;
}
}
free(ThisElement);
/* Read row index array: */
ThisElement = (char *) malloc(Indwidth+1);
if ( ThisElement == NULL ) IOHBTerminate("Insufficient memory for ThisElement.");
*(ThisElement+Indwidth) = (char) NULL;