-
Notifications
You must be signed in to change notification settings - Fork 130
/
sunmatrix_sparse.c
1236 lines (1002 loc) · 33.8 KB
/
sunmatrix_sparse.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
/*
* -----------------------------------------------------------------
* Programmer(s): Daniel Reynolds @ SMU
* David Gardner @ LLNL
* Based on code sundials_sparse.c by: Carol Woodward and
* Slaven Peles @ LLNL, and Daniel R. Reynolds @ SMU
* -----------------------------------------------------------------
* SUNDIALS Copyright Start
* Copyright (c) 2002-2022, Lawrence Livermore National Security
* and Southern Methodist University.
* All rights reserved.
*
* See the top-level LICENSE and NOTICE files for details.
*
* SPDX-License-Identifier: BSD-3-Clause
* SUNDIALS Copyright End
* -----------------------------------------------------------------
* This is the implementation file for the sparse implementation of
* the SUNMATRIX package.
* -----------------------------------------------------------------
*/
#include <stdio.h>
#include <stdlib.h>
#include <sunmatrix/sunmatrix_sparse.h>
#include <sundials/sundials_nvector.h>
#include <sundials/sundials_math.h>
#define ZERO RCONST(0.0)
#define ONE RCONST(1.0)
/* Private function prototypes */
static booleantype SMCompatible_Sparse(SUNMatrix A, SUNMatrix B);
static booleantype SMCompatible2_Sparse(SUNMatrix A, N_Vector x, N_Vector y);
static int Matvec_SparseCSC(SUNMatrix A, N_Vector x, N_Vector y);
static int Matvec_SparseCSR(SUNMatrix A, N_Vector x, N_Vector y);
static int format_convert(const SUNMatrix A, SUNMatrix B);
/*
* -----------------------------------------------------------------
* exported functions
* -----------------------------------------------------------------
*/
/*
* ==================================================================
* Private function prototypes (functions working on SlsMat)
* ==================================================================
*/
/* ----------------------------------------------------------------------------
* Function to create a new sparse matrix
*/
SUNMatrix SUNSparseMatrix(sunindextype M, sunindextype N,
sunindextype NNZ, int sparsetype,
SUNContext sunctx)
{
SUNMatrix A;
SUNMatrixContent_Sparse content;
/* return with NULL matrix on illegal input */
if ( (M <= 0) || (N <= 0) || (NNZ < 0) ) return(NULL);
if ( (sparsetype != CSC_MAT) && (sparsetype != CSR_MAT) ) return(NULL);
/* Create an empty matrix object */
A = NULL;
A = SUNMatNewEmpty(sunctx);
if (A == NULL) return(NULL);
/* Attach operations */
A->ops->getid = SUNMatGetID_Sparse;
A->ops->clone = SUNMatClone_Sparse;
A->ops->destroy = SUNMatDestroy_Sparse;
A->ops->zero = SUNMatZero_Sparse;
A->ops->copy = SUNMatCopy_Sparse;
A->ops->scaleadd = SUNMatScaleAdd_Sparse;
A->ops->scaleaddi = SUNMatScaleAddI_Sparse;
A->ops->matvec = SUNMatMatvec_Sparse;
A->ops->space = SUNMatSpace_Sparse;
/* Create content */
content = NULL;
content = (SUNMatrixContent_Sparse) malloc(sizeof *content);
if (content == NULL) { SUNMatDestroy(A); return(NULL); }
/* Attach content */
A->content = content;
/* Fill content */
content->sparsetype = sparsetype;
content->M = M;
content->N = N;
content->NNZ = NNZ;
switch(sparsetype){
case CSC_MAT:
content->NP = N;
content->rowvals = &(content->indexvals);
content->colptrs = &(content->indexptrs);
/* CSR indices */
content->colvals = NULL;
content->rowptrs = NULL;
break;
case CSR_MAT:
content->NP = M;
content->colvals = &(content->indexvals);
content->rowptrs = &(content->indexptrs);
/* CSC indices */
content->rowvals = NULL;
content->colptrs = NULL;
}
content->data = NULL;
content->indexvals = NULL;
content->indexptrs = NULL;
/* Allocate content */
content->data = (realtype *) calloc(NNZ, sizeof(realtype));
if (content->data == NULL) { SUNMatDestroy(A); return(NULL); }
content->indexvals = (sunindextype *) calloc(NNZ, sizeof(sunindextype));
if (content->indexvals == NULL) { SUNMatDestroy(A); return(NULL); }
content->indexptrs = (sunindextype *) calloc((content->NP + 1), sizeof(sunindextype));
if (content->indexptrs == NULL) { SUNMatDestroy(A); return(NULL); }
content->indexptrs[content->NP] = 0;
return(A);
}
/* ----------------------------------------------------------------------------
* Function to create a new sparse matrix from an existing dense matrix
* by copying all nonzero values into the sparse matrix structure. Returns NULL
* if the request for matrix storage cannot be satisfied.
*/
SUNMatrix SUNSparseFromDenseMatrix(SUNMatrix Ad, realtype droptol,
int sparsetype)
{
sunindextype i, j, nnz;
sunindextype M, N;
SUNMatrix As;
/* check for legal sparsetype, droptol and input matrix type */
if ( (sparsetype != CSR_MAT) && (sparsetype != CSC_MAT) )
return NULL;
if ( droptol < ZERO )
return NULL;
if (SUNMatGetID(Ad) != SUNMATRIX_DENSE)
return NULL;
/* set size of new matrix */
M = SM_ROWS_D(Ad);
N = SM_COLUMNS_D(Ad);
/* determine total number of nonzeros */
nnz = 0;
for (j=0; j<N; j++)
for (i=0; i<M; i++)
nnz += (SUNRabs(SM_ELEMENT_D(Ad,i,j)) > droptol);
/* allocate sparse matrix */
As = NULL;
As = SUNSparseMatrix(M, N, nnz, sparsetype, Ad->sunctx);
if (As == NULL) return NULL;
/* copy nonzeros from Ad into As, based on CSR/CSC type */
nnz = 0;
if (sparsetype == CSC_MAT) {
for (j=0; j<N; j++) {
(SM_INDEXPTRS_S(As))[j] = nnz;
for (i=0; i<M; i++) {
if ( SUNRabs(SM_ELEMENT_D(Ad,i,j)) > droptol ) {
(SM_INDEXVALS_S(As))[nnz] = i;
(SM_DATA_S(As))[nnz++] = SM_ELEMENT_D(Ad,i,j);
}
}
}
(SM_INDEXPTRS_S(As))[N] = nnz;
} else { /* CSR_MAT */
for (i=0; i<M; i++) {
(SM_INDEXPTRS_S(As))[i] = nnz;
for (j=0; j<N; j++) {
if ( SUNRabs(SM_ELEMENT_D(Ad,i,j)) > droptol ) {
(SM_INDEXVALS_S(As))[nnz] = j;
(SM_DATA_S(As))[nnz++] = SM_ELEMENT_D(Ad,i,j);
}
}
}
(SM_INDEXPTRS_S(As))[M] = nnz;
}
return(As);
}
/* ----------------------------------------------------------------------------
* Function to create a new sparse matrix from an existing band matrix
* by copying all nonzero values into the sparse matrix structure. Returns NULL
* if the request for matrix storage cannot be satisfied.
*/
SUNMatrix SUNSparseFromBandMatrix(SUNMatrix Ad, realtype droptol, int sparsetype)
{
sunindextype i, j, nnz;
sunindextype M, N;
SUNMatrix As;
/* check for legal sparsetype, droptol and input matrix type */
if ( (sparsetype != CSR_MAT) && (sparsetype != CSC_MAT) )
return NULL;
if ( droptol < ZERO )
return NULL;
if (SUNMatGetID(Ad) != SUNMATRIX_BAND)
return NULL;
/* set size of new matrix */
M = SM_ROWS_B(Ad);
N = SM_COLUMNS_B(Ad);
/* determine total number of nonzeros */
nnz = 0;
for (j=0; j<N; j++)
for (i=SUNMAX(0,j-SM_UBAND_B(Ad)); i<=SUNMIN(M-1,j+SM_LBAND_B(Ad)); i++)
nnz += (SUNRabs(SM_ELEMENT_B(Ad,i,j)) > droptol);
/* allocate sparse matrix */
As = SUNSparseMatrix(M, N, nnz, sparsetype, Ad->sunctx);
if (As == NULL) return NULL;
/* copy nonzeros from Ad into As, based on CSR/CSC type */
nnz = 0;
if (sparsetype == CSC_MAT) {
for (j=0; j<N; j++) {
(SM_INDEXPTRS_S(As))[j] = nnz;
for (i=SUNMAX(0,j-SM_UBAND_B(Ad)); i<=SUNMIN(M-1,j+SM_LBAND_B(Ad)); i++) {
if ( SUNRabs(SM_ELEMENT_B(Ad,i,j)) > droptol ) {
(SM_INDEXVALS_S(As))[nnz] = i;
(SM_DATA_S(As))[nnz++] = SM_ELEMENT_B(Ad,i,j);
}
}
}
(SM_INDEXPTRS_S(As))[N] = nnz;
} else { /* CSR_MAT */
for (i=0; i<M; i++) {
(SM_INDEXPTRS_S(As))[i] = nnz;
for (j=SUNMAX(0,i-SM_LBAND_B(Ad)); j<=SUNMIN(N-1,i+SM_UBAND_B(Ad)); j++) {
if ( SUNRabs(SM_ELEMENT_B(Ad,i,j)) > droptol ) {
(SM_INDEXVALS_S(As))[nnz] = j;
(SM_DATA_S(As))[nnz++] = SM_ELEMENT_B(Ad,i,j);
}
}
}
(SM_INDEXPTRS_S(As))[M] = nnz;
}
return(As);
}
/* ----------------------------------------------------------------------------
* Function to create a new CSR matrix from a CSC matrix.
*/
int SUNSparseMatrix_ToCSR(const SUNMatrix A, SUNMatrix* Bout)
{
if (A == NULL) return(SUNMAT_ILL_INPUT);
if (SM_SPARSETYPE_S(A) != CSC_MAT) return(SUNMAT_ILL_INPUT);
*Bout = SUNSparseMatrix(SM_ROWS_S(A), SM_COLUMNS_S(A), SM_NNZ_S(A), CSR_MAT, A->sunctx);
if (*Bout == NULL) return(SUNMAT_MEM_FAIL);
return format_convert(A, *Bout);
}
/* ----------------------------------------------------------------------------
* Function to create a new CSC matrix from a CSR matrix.
*/
int SUNSparseMatrix_ToCSC(const SUNMatrix A, SUNMatrix* Bout)
{
if (A == NULL) return(SUNMAT_ILL_INPUT);
if (SM_SPARSETYPE_S(A) != CSR_MAT) return(SUNMAT_ILL_INPUT);
*Bout = SUNSparseMatrix(SM_ROWS_S(A), SM_COLUMNS_S(A), SM_NNZ_S(A), CSC_MAT, A->sunctx);
if (*Bout == NULL) return(SUNMAT_MEM_FAIL);
return format_convert(A, *Bout);
}
/* ----------------------------------------------------------------------------
* Function to reallocate internal sparse matrix storage arrays so that the
* resulting sparse matrix holds indexptrs[NP] nonzeros. Returns 0 on success
* and 1 on failure (e.g. if A does not have sparse type, or if nnz is negative)
*/
int SUNSparseMatrix_Realloc(SUNMatrix A)
{
sunindextype nzmax;
/* check for valid matrix type */
if (SUNMatGetID(A) != SUNMATRIX_SPARSE)
return SUNMAT_ILL_INPUT;
/* get total number of nonzeros (return with failure if illegal) */
nzmax = (SM_INDEXPTRS_S(A))[SM_NP_S(A)];
if (nzmax < 0) return SUNMAT_ILL_INPUT;
/* perform reallocation */
SM_INDEXVALS_S(A) = (sunindextype *) realloc(SM_INDEXVALS_S(A), nzmax*sizeof(sunindextype));
SM_DATA_S(A) = (realtype *) realloc(SM_DATA_S(A), nzmax*sizeof(realtype));
SM_NNZ_S(A) = nzmax;
return SUNMAT_SUCCESS;
}
/* ----------------------------------------------------------------------------
* Function to reallocate internal sparse matrix storage arrays so that the
* resulting sparse matrix has storage for a specified number of nonzeros.
* Returns 0 on success and 1 on failure (e.g. if A does not have sparse type,
* or if nnz is negative)
*/
int SUNSparseMatrix_Reallocate(SUNMatrix A, sunindextype NNZ)
{
/* check for valid matrix type */
if (SUNMatGetID(A) != SUNMATRIX_SPARSE) return SUNMAT_ILL_INPUT;
/* check for valid nnz */
if (NNZ < 0) return SUNMAT_ILL_INPUT;
/* perform reallocation */
SM_INDEXVALS_S(A) = (sunindextype *) realloc(SM_INDEXVALS_S(A), NNZ*sizeof(sunindextype));
SM_DATA_S(A) = (realtype *) realloc(SM_DATA_S(A), NNZ*sizeof(realtype));
SM_NNZ_S(A) = NNZ;
return SUNMAT_SUCCESS;
}
/* ----------------------------------------------------------------------------
* Function to print the sparse matrix
*/
void SUNSparseMatrix_Print(SUNMatrix A, FILE* outfile)
{
sunindextype i, j;
char *matrixtype;
char *indexname;
/* should not be called unless A is a sparse matrix;
otherwise return immediately */
if (SUNMatGetID(A) != SUNMATRIX_SPARSE)
return;
/* perform operation */
if (SM_SPARSETYPE_S(A) == CSC_MAT) {
indexname = (char*) "col";
matrixtype = (char*) "CSC";
} else {
indexname = (char*) "row";
matrixtype = (char*) "CSR";
}
fprintf(outfile, "\n");
fprintf(outfile, "%ld by %ld %s matrix, NNZ: %ld \n",
(long int) SM_ROWS_S(A), (long int) SM_COLUMNS_S(A),
matrixtype, (long int) SM_NNZ_S(A));
for (j=0; j<SM_NP_S(A); j++) {
fprintf(outfile, "%s %ld : locations %ld to %ld\n", indexname,
(long int) j, (long int) (SM_INDEXPTRS_S(A))[j],
(long int) (SM_INDEXPTRS_S(A))[j+1]-1);
fprintf(outfile, " ");
for (i=(SM_INDEXPTRS_S(A))[j]; i<(SM_INDEXPTRS_S(A))[j+1]; i++) {
#if defined(SUNDIALS_EXTENDED_PRECISION)
fprintf(outfile, "%ld: %.32Lg ", (long int) (SM_INDEXVALS_S(A))[i],
(SM_DATA_S(A))[i]);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
fprintf(outfile, "%ld: %.16g ", (long int) (SM_INDEXVALS_S(A))[i],
(SM_DATA_S(A))[i]);
#else
fprintf(outfile, "%ld: %.8g ", (long int) (SM_INDEXVALS_S(A))[i],
(SM_DATA_S(A))[i]);
#endif
}
fprintf(outfile, "\n");
}
fprintf(outfile, "\n");
return;
}
/* ----------------------------------------------------------------------------
* Functions to access the contents of the sparse matrix structure
*/
sunindextype SUNSparseMatrix_Rows(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_ROWS_S(A);
else
return SUNMAT_ILL_INPUT;
}
sunindextype SUNSparseMatrix_Columns(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_COLUMNS_S(A);
else
return SUNMAT_ILL_INPUT;
}
sunindextype SUNSparseMatrix_NNZ(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_NNZ_S(A);
else
return SUNMAT_ILL_INPUT;
}
sunindextype SUNSparseMatrix_NP(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_NP_S(A);
else
return SUNMAT_ILL_INPUT;
}
int SUNSparseMatrix_SparseType(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_SPARSETYPE_S(A);
else
return SUNMAT_ILL_INPUT;
}
realtype* SUNSparseMatrix_Data(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_DATA_S(A);
else
return NULL;
}
sunindextype* SUNSparseMatrix_IndexValues(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_INDEXVALS_S(A);
else
return NULL;
}
sunindextype* SUNSparseMatrix_IndexPointers(SUNMatrix A)
{
if (SUNMatGetID(A) == SUNMATRIX_SPARSE)
return SM_INDEXPTRS_S(A);
else
return NULL;
}
/*
* -----------------------------------------------------------------
* implementation of matrix operations
* -----------------------------------------------------------------
*/
SUNMatrix_ID SUNMatGetID_Sparse(SUNMatrix A)
{
return SUNMATRIX_SPARSE;
}
SUNMatrix SUNMatClone_Sparse(SUNMatrix A)
{
SUNMatrix B = SUNSparseMatrix(SM_ROWS_S(A), SM_COLUMNS_S(A),
SM_NNZ_S(A), SM_SPARSETYPE_S(A), A->sunctx);
return(B);
}
void SUNMatDestroy_Sparse(SUNMatrix A)
{
if (A == NULL) return;
/* free content */
if (A->content != NULL) {
/* free data array */
if (SM_DATA_S(A)) {
free(SM_DATA_S(A));
SM_DATA_S(A) = NULL;
}
/* free index values array */
if (SM_INDEXVALS_S(A)) {
free(SM_INDEXVALS_S(A));
SM_INDEXVALS_S(A) = NULL;
SM_CONTENT_S(A)->rowvals = NULL;
SM_CONTENT_S(A)->colvals = NULL;
}
/* free index pointers array */
if (SM_INDEXPTRS_S(A)) {
free(SM_INDEXPTRS_S(A));
SM_INDEXPTRS_S(A) = NULL;
SM_CONTENT_S(A)->colptrs = NULL;
SM_CONTENT_S(A)->rowptrs = NULL;
}
/* free content struct */
free(A->content);
A->content = NULL;
}
/* free ops and matrix */
if (A->ops) { free(A->ops); A->ops = NULL; }
free(A); A = NULL;
return;
}
int SUNMatZero_Sparse(SUNMatrix A)
{
sunindextype i;
/* Perform operation */
for (i=0; i<SM_NNZ_S(A); i++) {
(SM_DATA_S(A))[i] = ZERO;
(SM_INDEXVALS_S(A))[i] = 0;
}
for (i=0; i<SM_NP_S(A); i++)
(SM_INDEXPTRS_S(A))[i] = 0;
(SM_INDEXPTRS_S(A))[SM_NP_S(A)] = 0;
return SUNMAT_SUCCESS;
}
int SUNMatCopy_Sparse(SUNMatrix A, SUNMatrix B)
{
sunindextype i, A_nz;
/* Verify that A and B are compatible */
if (!SMCompatible_Sparse(A, B))
return SUNMAT_ILL_INPUT;
/* Perform operation */
A_nz = (SM_INDEXPTRS_S(A))[SM_NP_S(A)];
/* ensure that B is allocated with at least as
much memory as we have nonzeros in A */
if (SM_NNZ_S(B) < A_nz) {
SM_INDEXVALS_S(B) = (sunindextype *) realloc(SM_INDEXVALS_S(B), A_nz*sizeof(sunindextype));
SM_DATA_S(B) = (realtype *) realloc(SM_DATA_S(B), A_nz*sizeof(realtype));
SM_NNZ_S(B) = A_nz;
}
/* zero out B so that copy works correctly */
if (SUNMatZero_Sparse(B) != SUNMAT_SUCCESS)
return SUNMAT_OPERATION_FAIL;
/* copy the data and row indices over */
for (i=0; i<A_nz; i++){
(SM_DATA_S(B))[i] = (SM_DATA_S(A))[i];
(SM_INDEXVALS_S(B))[i] = (SM_INDEXVALS_S(A))[i];
}
/* copy the column pointers over */
for (i=0; i<SM_NP_S(A); i++) {
(SM_INDEXPTRS_S(B))[i] = (SM_INDEXPTRS_S(A))[i];
}
(SM_INDEXPTRS_S(B))[SM_NP_S(A)] = A_nz;
return SUNMAT_SUCCESS;
}
int SUNMatScaleAddI_Sparse(realtype c, SUNMatrix A)
{
sunindextype j, i, p, nz, newvals, M, N, cend;
booleantype newmat, found;
sunindextype *w, *Ap, *Ai, *Cp, *Ci;
realtype *x, *Ax, *Cx;
SUNMatrix C;
/* store shortcuts to matrix dimensions (M is inner dimension, N is outer) */
if (SM_SPARSETYPE_S(A) == CSC_MAT) {
M = SM_ROWS_S(A);
N = SM_COLUMNS_S(A);
}
else {
M = SM_COLUMNS_S(A);
N = SM_ROWS_S(A);
}
/* access data arrays from A (return if failure) */
Ap = Ai = NULL;
Ax = NULL;
if (SM_INDEXPTRS_S(A)) Ap = SM_INDEXPTRS_S(A);
else return (SUNMAT_MEM_FAIL);
if (SM_INDEXVALS_S(A)) Ai = SM_INDEXVALS_S(A);
else return (SUNMAT_MEM_FAIL);
if (SM_DATA_S(A)) Ax = SM_DATA_S(A);
else return (SUNMAT_MEM_FAIL);
/* determine if A: contains values on the diagonal (so I can just be added in);
if not, then increment counter for extra storage that should be required. */
newvals = 0;
for (j=0; j < SUNMIN(M,N); j++) {
/* scan column (row if CSR) of A, searching for diagonal value */
found = SUNFALSE;
for (i=Ap[j]; i<Ap[j+1]; i++) {
if (Ai[i] == j) {
found = SUNTRUE;
break;
}
}
/* if no diagonal found, increment necessary storage counter */
if (!found) newvals += 1;
}
/* If extra nonzeros required, check whether matrix has sufficient storage space
for new nonzero entries (so I can be inserted into existing storage) */
newmat = SUNFALSE; /* no reallocation needed */
if (newvals > (SM_NNZ_S(A) - Ap[N]))
newmat = SUNTRUE;
/* perform operation based on existing/necessary structure */
/* case 1: A already contains a diagonal */
if (newvals == 0) {
/* iterate through columns, adding 1.0 to diagonal */
for (j=0; j < SUNMIN(M,N); j++)
for (i=Ap[j]; i<Ap[j+1]; i++)
if (Ai[i] == j) {
Ax[i] = ONE + c*Ax[i];
} else {
Ax[i] = c*Ax[i];
}
/* case 2: A has sufficient storage, but does not already contain a diagonal */
} else if (!newmat) {
/* create work arrays for nonzero indices and values in a single column (row) */
w = (sunindextype *) malloc(M * sizeof(sunindextype));
x = (realtype *) malloc(M * sizeof(realtype));
/* determine storage location where last column (row) should end */
nz = Ap[N] + newvals;
/* store pointer past last column (row) from original A,
and store updated value in revised A */
cend = Ap[N];
Ap[N] = nz;
/* iterate through columns (rows) backwards */
for (j=N-1; j>=0; j--) {
/* clear out temporary arrays for this column (row) */
for (i=0; i<M; i++) {
w[i] = 0;
x[i] = RCONST(0.0);
}
/* iterate down column (row) of A, collecting nonzeros */
for (p=Ap[j]; p<cend; p++) {
w[Ai[p]] += 1; /* indicate that row (column) is filled */
x[Ai[p]] = c*Ax[p]; /* collect/scale value */
}
/* add identity to this column (row) */
if (j < M) {
w[j] += 1; /* indicate that row (column) is filled */
x[j] += ONE; /* update value */
}
/* fill entries of A with this column's (row's) data */
for (i=M-1; i>=0; i--) {
if ( w[i] > 0 ) {
Ai[--nz] = i;
Ax[nz] = x[i];
}
}
/* store ptr past this col (row) from orig A, update value for new A */
cend = Ap[j];
Ap[j] = nz;
}
/* clean up */
free(w);
free(x);
/* case 3: A must be reallocated with sufficient storage */
} else {
/* create work arrays for nonzero indices and values */
w = (sunindextype *) malloc(M * sizeof(sunindextype));
x = (realtype *) malloc(M * sizeof(realtype));
/* create new matrix for sum */
C = SUNSparseMatrix(SM_ROWS_S(A), SM_COLUMNS_S(A),
Ap[N] + newvals,
SM_SPARSETYPE_S(A), A->sunctx);
/* access data from CSR structures (return if failure) */
Cp = Ci = NULL;
Cx = NULL;
if (SM_INDEXPTRS_S(C)) Cp = SM_INDEXPTRS_S(C);
else return (SUNMAT_MEM_FAIL);
if (SM_INDEXVALS_S(C)) Ci = SM_INDEXVALS_S(C);
else return (SUNMAT_MEM_FAIL);
if (SM_DATA_S(C)) Cx = SM_DATA_S(C);
else return (SUNMAT_MEM_FAIL);
/* initialize total nonzero count */
nz = 0;
/* iterate through columns (rows for CSR) */
for (j=0; j<N; j++) {
/* set current column (row) pointer to current # nonzeros */
Cp[j] = nz;
/* clear out temporary arrays for this column (row) */
for (i=0; i<M; i++) {
w[i] = 0;
x[i] = 0.0;
}
/* iterate down column (along row) of A, collecting nonzeros */
for (p=Ap[j]; p<Ap[j+1]; p++) {
w[Ai[p]] += 1; /* indicate that row is filled */
x[Ai[p]] = c*Ax[p]; /* collect/scale value */
}
/* add identity to this column (row) */
if (j < M) {
w[j] += 1; /* indicate that row is filled */
x[j] += ONE; /* update value */
}
/* fill entries of C with this column's (row's) data */
for (i=0; i<M; i++) {
if ( w[i] > 0 ) {
Ci[nz] = i;
Cx[nz++] = x[i];
}
}
}
/* indicate end of data */
Cp[N] = nz;
/* update A's structure with C's values; nullify C's pointers */
SM_NNZ_S(A) = SM_NNZ_S(C);
if (SM_DATA_S(A))
free(SM_DATA_S(A));
SM_DATA_S(A) = SM_DATA_S(C);
SM_DATA_S(C) = NULL;
if (SM_INDEXVALS_S(A))
free(SM_INDEXVALS_S(A));
SM_INDEXVALS_S(A) = SM_INDEXVALS_S(C);
SM_INDEXVALS_S(C) = NULL;
if (SM_INDEXPTRS_S(A))
free(SM_INDEXPTRS_S(A));
SM_INDEXPTRS_S(A) = SM_INDEXPTRS_S(C);
SM_INDEXPTRS_S(C) = NULL;
/* clean up */
SUNMatDestroy_Sparse(C);
free(w);
free(x);
}
return SUNMAT_SUCCESS;
}
int SUNMatScaleAdd_Sparse(realtype c, SUNMatrix A, SUNMatrix B)
{
sunindextype j, i, p, nz, newvals, M, N, cend;
booleantype newmat;
sunindextype *w, *Ap, *Ai, *Bp, *Bi, *Cp, *Ci;
realtype *x, *Ax, *Bx, *Cx;
SUNMatrix C;
/* Verify that A and B are compatible */
if (!SMCompatible_Sparse(A, B))
return SUNMAT_ILL_INPUT;
/* store shortcuts to matrix dimensions (M is inner dimension, N is outer) */
if (SM_SPARSETYPE_S(A) == CSC_MAT) {
M = SM_ROWS_S(A);
N = SM_COLUMNS_S(A);
}
else {
M = SM_COLUMNS_S(A);
N = SM_ROWS_S(A);
}
/* access data arrays from A and B (return if failure) */
Ap = Ai = Bp = Bi = NULL;
Ax = Bx = NULL;
if (SM_INDEXPTRS_S(A)) Ap = SM_INDEXPTRS_S(A);
else return(SUNMAT_MEM_FAIL);
if (SM_INDEXVALS_S(A)) Ai = SM_INDEXVALS_S(A);
else return(SUNMAT_MEM_FAIL);
if (SM_DATA_S(A)) Ax = SM_DATA_S(A);
else return(SUNMAT_MEM_FAIL);
if (SM_INDEXPTRS_S(B)) Bp = SM_INDEXPTRS_S(B);
else return(SUNMAT_MEM_FAIL);
if (SM_INDEXVALS_S(B)) Bi = SM_INDEXVALS_S(B);
else return(SUNMAT_MEM_FAIL);
if (SM_DATA_S(B)) Bx = SM_DATA_S(B);
else return(SUNMAT_MEM_FAIL);
/* create work arrays for row indices and nonzero column values */
w = (sunindextype *) malloc(M * sizeof(sunindextype));
x = (realtype *) malloc(M * sizeof(realtype));
/* determine if A already contains the sparsity pattern of B */
newvals = 0;
for (j=0; j<N; j++) {
/* clear work array */
for (i=0; i<M; i++) w[i] = 0;
/* scan column of A, incrementing w by one */
for (i=Ap[j]; i<Ap[j+1]; i++)
w[Ai[i]] += 1;
/* scan column of B, decrementing w by one */
for (i=Bp[j]; i<Bp[j+1]; i++)
w[Bi[i]] -= 1;
/* if any entry of w is negative, A doesn't contain B's sparsity,
so increment necessary storage counter */
for (i=0; i<M; i++)
if (w[i] < 0) newvals += 1;
}
/* If extra nonzeros required, check whether A has sufficient storage space
for new nonzero entries (so B can be inserted into existing storage) */
newmat = SUNFALSE; /* no reallocation needed */
if (newvals > (SM_NNZ_S(A) - Ap[N]))
newmat = SUNTRUE;
/* perform operation based on existing/necessary structure */
/* case 1: A already contains sparsity pattern of B */
if (newvals == 0) {
/* iterate through columns, adding matrices */
for (j=0; j<N; j++) {
/* clear work array */
for (i=0; i<M; i++)
x[i] = ZERO;
/* scan column of B, updating work array */
for (i = Bp[j]; i < Bp[j+1]; i++)
x[Bi[i]] = Bx[i];
/* scan column of A, updating array entries appropriately */
for (i = Ap[j]; i < Ap[j+1]; i++)
Ax[i] = c*Ax[i] + x[Ai[i]];
}
/* case 2: A has sufficient storage, but does not already contain B's sparsity */
} else if (!newmat) {
/* determine storage location where last column (row) should end */
nz = Ap[N] + newvals;
/* store pointer past last column (row) from original A,
and store updated value in revised A */
cend = Ap[N];
Ap[N] = nz;
/* iterate through columns (rows) backwards */
for (j=N-1; j>=0; j--) {
/* clear out temporary arrays for this column (row) */
for (i=0; i<M; i++) {
w[i] = 0;
x[i] = RCONST(0.0);
}
/* iterate down column (row) of A, collecting nonzeros */
for (p=Ap[j]; p<cend; p++) {
w[Ai[p]] += 1; /* indicate that row (column) is filled */
x[Ai[p]] = c*Ax[p]; /* collect/scale value */
}
/* iterate down column of B, collecting nonzeros */
for (p=Bp[j]; p<Bp[j+1]; p++) {
w[Bi[p]] += 1; /* indicate that row is filled */
x[Bi[p]] += Bx[p]; /* collect value */
}
/* fill entries of A with this column's (row's) data */
for (i=M-1; i>=0; i--) {
if ( w[i] > 0 ) {
Ai[--nz] = i;
Ax[nz] = x[i];
}
}
/* store ptr past this col (row) from orig A, update value for new A */
cend = Ap[j];
Ap[j] = nz;
}
/* case 3: A must be reallocated with sufficient storage */
} else {
/* create new matrix for sum */
C = SUNSparseMatrix(SM_ROWS_S(A), SM_COLUMNS_S(A),
Ap[N] + newvals, SM_SPARSETYPE_S(A), A->sunctx);
/* access data from CSR structures (return if failure) */
Cp = Ci = NULL;
Cx = NULL;
if (SM_INDEXPTRS_S(C)) Cp = SM_INDEXPTRS_S(C);
else return(SUNMAT_MEM_FAIL);
if (SM_INDEXVALS_S(C)) Ci = SM_INDEXVALS_S(C);
else return(SUNMAT_MEM_FAIL);
if (SM_DATA_S(C)) Cx = SM_DATA_S(C);
else return(SUNMAT_MEM_FAIL);
/* initialize total nonzero count */
nz = 0;
/* iterate through columns (rows) */
for (j=0; j<N; j++) {
/* set current column (row) pointer to current # nonzeros */
Cp[j] = nz;
/* clear out temporary arrays for this column (row) */
for (i=0; i<M; i++) {
w[i] = 0;
x[i] = RCONST(0.0);
}
/* iterate down column of A, collecting nonzeros */
for (p=Ap[j]; p<Ap[j+1]; p++) {
w[Ai[p]] += 1; /* indicate that row is filled */
x[Ai[p]] = c*Ax[p]; /* collect/scale value */
}
/* iterate down column of B, collecting nonzeros */
for (p=Bp[j]; p<Bp[j+1]; p++) {
w[Bi[p]] += 1; /* indicate that row is filled */
x[Bi[p]] += Bx[p]; /* collect value */
}
/* fill entries of C with this column's data */
for (i=0; i<M; i++) {
if ( w[i] > 0 ) {
Ci[nz] = i;
Cx[nz++] = x[i];
}
}
}
/* indicate end of data */
Cp[N] = nz;
/* update A's structure with C's values; nullify C's pointers */
SM_NNZ_S(A) = SM_NNZ_S(C);
free(SM_DATA_S(A));
SM_DATA_S(A) = SM_DATA_S(C);
SM_DATA_S(C) = NULL;
free(SM_INDEXVALS_S(A));
SM_INDEXVALS_S(A) = SM_INDEXVALS_S(C);
SM_INDEXVALS_S(C) = NULL;
free(SM_INDEXPTRS_S(A));
SM_INDEXPTRS_S(A) = SM_INDEXPTRS_S(C);
SM_INDEXPTRS_S(C) = NULL;
/* clean up */
SUNMatDestroy_Sparse(C);
}