-
Notifications
You must be signed in to change notification settings - Fork 0
/
levmar.c
2242 lines (1987 loc) · 68.7 KB
/
levmar.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
/////////////////////////////////////////////////////////////////////////////////
//
// Sparse levenberg marquardt algorithm. The algorithm is similar to
// minpack's lmdif algorithm except that optimizations for sparse matrices
// are used.
//
// Copyright (C) 2021 Florian Königstein
//
// The code was developed by taking minpack's lmdif algorithm and replacing
// where necessary the dense matrix operations by sparse matrix operations.
// The algorithm depends on the SuiteSparse project
// (by Christopher Lourenco, JinHao Chen, Erick Moreno-Centeno,
// and Timothy A.Davis) and optionally on the NVIDIA GPU Computing Toolkit
// (copyright NVIDIA, 2701 San TomasExpressway, Santa Clara, CA 95050).
// Some code and ideas are borrowed from the the sparseLM library
// (by Manolis Lourakis).
// Copyright notice of minpack: march 1980 argonne national laboratory.
// Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More
//
// This program is free software; you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation; either version 2 of the License, or
// (at your option) any later version.
//
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU General Public License for more details.
//
/////////////////////////////////////////////////////////////////////////////////
//#define WITH_NVIDIA_CUDA_GPU_SUPPORT
#include "levmar.h"
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <memory.h>
#include <inttypes.h>
#include <cholmod.h>
#include <SuiteSparseQR_C.h>
#ifdef WITH_NVIDIA_CUDA_GPU_SUPPORT
#include <cuda_runtime.h>
#include "cuda.h"
#endif
static const double MACHEP = 1.2e-16;
/* smallest nonzero number */
static const double DWARF = 1.0e-38; // in enorm 3.834e-20
static const double RGIANT = 1.304e19;
#ifdef WITH_DEBUG_PRINT_MATRIX_VECTOR
void print_matrix_C(char* str, double* M, int64_t n, int64_t m)
{
int64_t i, j;
printf("%s :\n", str);
for (i = 0; i < n; i++)
{
for (j = 0; j < m; j++)
printf("%.16f ", M[i + n * j]);
printf("\n");
}
printf("\n");
}
void print_vector_C(char* str, double* v, int64_t m)
{
int64_t j;
printf("%s :\n", str);
for (j = 0; j < m; j++)
printf("%.16f ", v[j]);
printf("\n\n");
}
void print_int_vector_C(char* str, int64_t* v, int64_t m)
{
int64_t j;
printf("%s :\n", str);
for (j = 0; j < m; j++)
printf("%lld ", v[j]);
printf("\n\n");
}
#endif
#ifdef WITH_NVIDIA_CUDA_GPU_SUPPORT
int cudaInitDevice(CUcontext* cuContext)
{
int selectedDevice = 0;
int deviceCount;
cudaGetDeviceCount(&deviceCount);
if (deviceCount == 0)
{
return 0;
}
if (selectedDevice >= deviceCount)
{
return 0;
}
CUdevice cuDevice;
return CUDA_SUCCESS == cuDeviceGet(&cuDevice, selectedDevice) &&
CUDA_SUCCESS == cuCtxCreate(cuContext, CU_CTX_MAP_HOST | CU_CTX_BLOCKING_SYNC, cuDevice);
}
void cudaReleaseDevice(CUcontext* cuContext)
{
CUresult cerr;
cerr = cuCtxDestroy(*cuContext);
}
#endif
static double enorm(int64_t n, double x[])
{
/*
* **********
*
* function enorm
*
* given an n-vector x, this function calculates the
* euclidean norm of x.
*
* the euclidean norm is computed by accumulating the sum of
* squares in three different sums. the sums of squares for the
* small and large components are scaled so that no overflows
* occur. non-destructive underflows are permitted. underflows
* and overflows do not occur in the computation of the unscaled
* sum of squares for the intermediate components.
* the definitions of small, intermediate and large components
* depend on two constants, rdwarf and rgiant. the main
* restrictions on these constants are that rdwarf**2 not
* underflow and rgiant**2 not overflow. the constants
* given here are suitable for every known computer.
*
* the function statement is
*
* double precision function enorm(n,x)
*
* where
*
* n is a positive integer input variable.
*
* x is an input array of length n.
*
* subprograms called
*
* fortran-supplied ... dabs,dsqrt
*
* argonne national laboratory. minpack project. march 1980.
* burton s. garbow, kenneth e. hillstrom, jorge j. more
*
* **********
*/
int64_t i;
double agiant, s1, s2, s3, xabs, x1max, x3max;
double temp;
s1 = 0.0;
s2 = 0.0;
s3 = 0.0;
x1max = 0.0;
x3max = 0.0;
agiant = RGIANT / ((double)n);
for (i = 0; i < n; i++)
{
xabs = fabs(x[i]);
if ((xabs > DWARF) && (xabs <agiant))
{
/*
* sum for intermediate components.
*/
s2 += xabs * xabs;
continue;
}
if (xabs > DWARF)
{
/*
* sum for large components.
*/
if (xabs > x1max)
{
temp = x1max / xabs;
s1 = 1.0 + s1 * temp * temp;
x1max = xabs;
}
else
{
temp = xabs / x1max;
s1 += temp * temp;
}
continue;
}
/*
* sum for small components.
*/
if (xabs > x3max)
{
temp = x3max / xabs;
s3 = 1.0 + s3 * temp * temp;
x3max = xabs;
}
else
{
if (xabs != 0.0)
{
temp = xabs / x3max;
s3 += temp * temp;
}
}
}
/*
* calculation of norm.
*/
if (s1 != 0.0)
{
return x1max * sqrt(s1 + (s2 / x1max) / x1max);
}
if (s2 != 0.0)
{
if (s2 >= x3max)
temp = s2 * (1.0 + (x3max / s2) * (x3max * s3));
else
temp = x3max * ((s2 / x3max) + (x3max * s3));
return sqrt(temp);
}
else
{
return x3max * sqrt(s3);
}
}
/* convert a matrix from the CRS format to CCS. If crs->val is NULL, only the nonzero pattern is converted */
/* ccs must already have been allocated and have the same nr, nc and nnz */
SuiteSparse_long splm_crsm2ccsm(splm_crsm* crs, splm_ccsm* ccs)
{
register SuiteSparse_long i, j, k, l;
SuiteSparse_long nr, nc, nnz, jmax;
SuiteSparse_long* colidx, * rowptr, * rowidx, * colptr;
SuiteSparse_long* colcounts; // counters for the number of nonzeros in each column
nr = crs->nr; nc = crs->nc;
nnz = crs->nnz;
if (0 == (colcounts = (SuiteSparse_long*)malloc(nc * sizeof(SuiteSparse_long))))
{
return -1;
}
ccs->nr = nr; ccs->nc = nc; // ensure that ccs has the correct dimensions
colidx = crs->colidx; rowptr = crs->rowptr;
rowidx = ccs->rowidx; colptr = ccs->colptr;
for (j = 0; j < nc; j++)
colcounts[j] = 0;
/* 1st pass: count #nonzeros in each column */
for (j = rowptr[nr]; 0 != j; )
++(colcounts[colidx[--j]]);
/* 2nd pass: copy every nonzero to its right position into the CCS structure */
for (j = k = 0; j < nc; ++j)
{
colptr[j] = k;
k += colcounts[j];
colcounts[j] = 0; // clear to avoid separate loop below
}
colptr[nc] = nnz;
/* colcounts[j] will count the #nonzeros in col. j seen before the current row */
if (crs->val) { // are there any values to copy?
register double* crsv, * ccsv;
crsv = crs->val; ccsv = ccs->val;
for (i = 0; i < nr; ++i)
{
jmax = rowptr[i + 1];
for (j = rowptr[i]; j < jmax; ++j)
{
l = colidx[j];
k = colptr[l];
k += colcounts[l]++;
rowidx[k] = i;
/* copy values */
ccsv[k] = crsv[j];
}
}
}
else
{ // no values, copy just structure
for (i = 0; i < nr; ++i)
{
jmax = rowptr[i + 1];
for (j = rowptr[i]; j < jmax; ++j)
{
l = colidx[j];
k = colptr[l];
k += colcounts[l]++;
rowidx[k] = i;
}
}
}
free(colcounts);
return 0;
}
void splm_ccsm_init_invalid(splm_ccsm* sm)
{
sm->val = 0;
sm->rowidx = 0;
sm->colptr = 0;
sm->nr = sm->nc = sm->nnz = -1;
sm->cm_owner = 0;
sm->cm_common = 0;
}
/* allocate a sparse CCS matrix */
SuiteSparse_long splm_ccsm_alloc(splm_ccsm* sm, SuiteSparse_long nr, SuiteSparse_long nc, SuiteSparse_long nnz)
{
sm->val = (double*)malloc(nnz * sizeof(double));
sm->rowidx = (SuiteSparse_long*)malloc(nnz * sizeof(SuiteSparse_long));
sm->colptr = (SuiteSparse_long*)malloc((nc + 1) * sizeof(SuiteSparse_long));
if (!sm->val || !sm->rowidx || !sm->colptr)
{
if (sm->val) { free(sm->val); sm->val = 0; }
if (sm->rowidx) { free(sm->rowidx); sm->rowidx = 0; }
if (sm->colptr) { free(sm->colptr); sm->colptr = 0; }
sm->nr = sm->nc = sm->nnz = -1;
return -1;
}
sm->nr = nr;
sm->nc = nc;
sm->nnz = nnz;
sm->cm_owner = 0;
sm->cm_common = 0;
return 0;
}
/* free a sparse CCS matrix (but NOT the cholmod part if any)*/
void splm_ccsm_free(splm_ccsm* sm)
{
if (sm->val) { free(sm->val); sm->val = 0; }
if (sm->rowidx) { free(sm->rowidx); sm->rowidx = 0; }
if (sm->colptr) { free(sm->colptr); sm->colptr = 0; }
sm->nr = sm->nc = sm->nnz = -1;
}
/* returns the index of the (i, j) element. No bounds checking! */
SuiteSparse_long splm_ccsm_elmidx(splm_ccsm* sm, SuiteSparse_long i, SuiteSparse_long j)
{
register SuiteSparse_long low, high, mid, * Rowidx, diff;
low = sm->colptr[j];
high = sm->colptr[j + 1] - 1;
Rowidx = sm->rowidx;
/* binary search for finding the element at row i */
while (low <= high) {
//if(i<Rowidx[low] || i>Rowidx[high]) return -1; /* not found */
mid = (low + high) >> 1; //(low+high)/2;
//mid=low+((high-low)>>1); /* ensures no index overflows */
diff = i - Rowidx[mid];
if (diff < 0)
high = mid - 1;
else if (diff > 0)
low = mid + 1;
else
return mid;
}
return -1; /* not found */
}
/* returns the maximum number of nonzero elements across all columns */
SuiteSparse_long splm_ccsm_col_maxnelms(splm_ccsm* sm)
{
register SuiteSparse_long j, n, max;
for (j = sm->nc, max = -1; j-- > 0; )
if ((n = sm->colptr[j + 1] - sm->colptr[j]) > max) max = n;
return max;
}
/* returns the number of nonzero elements in col j and
* fills up the vidxs and iidxs arrays with the val and row
* indexes of the elements found, respectively.
* vidxs and iidxs are assumed preallocated and of max. size sm->nr
*/
SuiteSparse_long splm_ccsm_col_elmidxs(splm_ccsm* sm, SuiteSparse_long j, SuiteSparse_long* vidxs, SuiteSparse_long* iidxs)
{
register SuiteSparse_long i, k;
SuiteSparse_long low, high, * rowidx = sm->rowidx;
low = sm->colptr[j];
high = sm->colptr[j + 1];
for (i = low, k = 0; i < high; ++i, ++k) {
vidxs[k] = i;
iidxs[k] = rowidx[i];
}
return k;
}
void splm_crsm_init_invalid(splm_crsm* sm)
{
sm->val = 0;
sm->colidx = 0;
sm->rowptr = 0;
sm->nr = sm->nc = sm->nnz = -1;
}
/* allocate all fields except values */
SuiteSparse_long splm_crsm_alloc_novalues(splm_crsm* sm, SuiteSparse_long nr, SuiteSparse_long nc, SuiteSparse_long nnz)
{
sm->nr = nr;
sm->nc = nc;
sm->nnz = nnz;
sm->val = 0;
sm->colidx = nnz > 0 ? (SuiteSparse_long*)malloc(nnz * sizeof(SuiteSparse_long)) : 0;
sm->rowptr = (SuiteSparse_long*)malloc((nr + 1) * sizeof(SuiteSparse_long));
if ((nnz > 0 && !sm->colidx) || !sm->rowptr) {
if(sm->colidx) { free(sm->colidx); sm->colidx=0; }
if (sm->rowptr) { free(sm->rowptr); sm->rowptr = 0; }
sm->nr = sm->nc = sm->nnz = -1;
return -1;
}
return 0;
}
/* allocate val and colidx after allocating rowptr with splm_crsm_alloc_novalues(...,0) */
SuiteSparse_long splm_crsm_alloc_rest(splm_crsm* sm, SuiteSparse_long nnz)
{
if(sm->nr<0 || sm->nc<0 || 0==sm->rowptr)
return -1;
sm->nnz = nnz;
sm->val = (double*)malloc(nnz * sizeof(double));;
sm->colidx = (SuiteSparse_long*)malloc(nnz * sizeof(SuiteSparse_long));
if (!sm->val || !sm->colidx)
{
if (sm->val) { free(sm->val); sm->val = 0; }
if (sm->colidx) { free(sm->colidx); sm->colidx = 0; }
if (sm->rowptr) { free(sm->rowptr); sm->rowptr = 0; }
sm->nr = sm->nc = sm->nnz = -1;
return -1;
}
return 0;
}
/* free a sparse CRS matrix */
void splm_crsm_free(splm_crsm* sm)
{
if (sm->val) { free(sm->val); sm->val = 0; }
if (sm->colidx) { free(sm->colidx); sm->colidx = 0; }
if (sm->rowptr) { free(sm->rowptr); sm->rowptr = 0; }
sm->nr = sm->nc = sm->nnz = -1;
}
splm_ccsm* cholmod_sparse_to_splm_ccsm(cholmod_sparse* cmsp, cholmod_common* cm_common)
{ /* If this function succeeds (returns != 0), you must not call cholmod_free_sparse(cmsp) on destroying.
Instead you must call splm_ccsm_destruct(ccsm) where ccsm is the pointer returned by
cholmod_sparse_to_splm_ccsm). */
splm_ccsm* sm;
if (!cmsp || !cmsp->packed || !cmsp->sorted || cmsp->dtype != CHOLMOD_DOUBLE || cmsp->xtype != CHOLMOD_REAL)
{
return 0;
}
if (0 != (sm = malloc(sizeof(splm_ccsm))))
{
sm->nr = cmsp->nrow;
sm->nc = cmsp->ncol;
sm->nnz = cmsp->nzmax;
sm->colptr = (SuiteSparse_long*)cmsp->p;
sm->rowidx = (SuiteSparse_long*)cmsp->i;
sm->val = (double*)cmsp->x;
sm->cm_owner = cmsp;
sm->cm_common = cm_common;
}
return sm;
}
void splm_ccsm_destruct(splm_ccsm* sm)
{
if (sm)
{
if (sm->cm_owner)
{
cholmod_l_free_sparse(&sm->cm_owner, sm->cm_common);
}
else
{
splm_ccsm_free(sm);
}
free(sm);
}
}
cholmod_sparse* cholmod_sparse_mirror_from_ccsm(splm_ccsm* sm)
{ /* You must not call cholmod_free_sparse() for the maxtix returned by this function. */
if (!sm)
{
return 0;
}
sm->static_cm_mirror.nrow = sm->nr;
sm->static_cm_mirror.ncol = sm->nc;
sm->static_cm_mirror.nzmax = sm->nnz;
sm->static_cm_mirror.x = sm->val;
sm->static_cm_mirror.p = sm->colptr;
sm->static_cm_mirror.i = sm->rowidx;
sm->static_cm_mirror.nz = sm->static_cm_mirror.z = 0;
sm->static_cm_mirror.stype = 0; // unsymmetric matrix
sm->static_cm_mirror.itype = CHOLMOD_LONG;
sm->static_cm_mirror.xtype = CHOLMOD_REAL;
sm->static_cm_mirror.dtype = CHOLMOD_DOUBLE;
sm->static_cm_mirror.sorted = sm->static_cm_mirror.packed = 1;
return &(sm->static_cm_mirror);
}
SuiteSparse_long splm_SuiteSparseQR
(
// inputs, not modified
int ordering, // all, except 3:given treated as 0:fixed
double tol, // columns with 2-norm <= tol are treated as zero
SuiteSparse_long econ, // number of rows of C and R to return; a value
// less than the rank r of A is treated as r, and
// a value greater than m is treated as m.
// That is, e = max(min(m,econ),rank(A)) gives the
// number of rows of C and R, and columns of C'.
int getCTX, // if 0: return Z = C of size econ-by-bncols
// if 1: return Z = C' of size bncols-by-econ
// if 2: return Z = X of size econ-by-bncols
splm_ccsm* A, // m-by-n sparse matrix
// B is either sparse or dense. If Bsparse is non-NULL, B is sparse and
// Bdense is ignored. If Bsparse is NULL and Bdense is non-NULL, then B is
// dense. B is not present if both are NULL.
splm_ccsm* Bsparse, // B is m-by-bncols
cholmod_dense* Bdense,
// output arrays, neither allocated nor defined on input.
// Z is the matrix C, C', or X, either sparse or dense. If p_Zsparse is
// non-NULL, then Z is returned as a sparse matrix. If p_Zsparse is NULL
// and p_Zdense is non-NULL, then Z is returned as a dense matrix. Neither
// are returned if both arguments are NULL.
splm_ccsm** p_Zsparse,
cholmod_dense** p_Zdense,
splm_ccsm** p_R, // the R factor
SuiteSparse_long** p_E, // size n; fill-reducing ordering of A.
splm_ccsm** p_H, // the Householder vectors (m-by-nh)
SuiteSparse_long** p_HPinv, // size m; row permutation for H
cholmod_dense** p_HTau, // size 1-by-nh, Householder coefficients
// workspace and parameters
cholmod_common* cc
)
{
cholmod_sparse* p_cmZsparse = 0;
cholmod_sparse* p_cmR = 0;
cholmod_sparse* p_cmH = 0;
SuiteSparse_long rank = SuiteSparseQR_C(ordering, tol, econ, getCTX,
cholmod_sparse_mirror_from_ccsm(A),
cholmod_sparse_mirror_from_ccsm(Bsparse), Bdense,
p_Zsparse ? (&p_cmZsparse) : 0, p_Zdense,
p_R ? (&p_cmR) : 0, p_E,
p_H ? (&p_cmH) : 0, p_HPinv, p_HTau, cc);
if(cc->useGPU && CHOLMOD_GPU_PROBLEM == cc->status)
{
cc->useGPU = 0;
rank = SuiteSparseQR_C(ordering, tol, econ, getCTX,
cholmod_sparse_mirror_from_ccsm(A),
cholmod_sparse_mirror_from_ccsm(Bsparse), Bdense,
p_Zsparse ? (&p_cmZsparse) : 0, p_Zdense,
p_R ? (&p_cmR) : 0, p_E,
p_H ? (&p_cmH) : 0, p_HPinv, p_HTau, cc);
}
if (cc->status < CHOLMOD_OK ||
(p_Zsparse && !p_cmZsparse) || (p_R && !(p_cmR)) || (p_H && !(p_cmH)))
{
rank = -1;
}
else if (p_cmZsparse && 0 == (*p_Zsparse = cholmod_sparse_to_splm_ccsm(p_cmZsparse, cc)))
{
rank = -1;
}
else if (p_cmR && 0 == (*p_R = cholmod_sparse_to_splm_ccsm(p_cmR, cc)))
{
rank = -1;
}
else if (p_H && 0 == (*p_H = cholmod_sparse_to_splm_ccsm(p_cmH, cc)))
{
rank = -1;
}
if (rank < 0)
{
if (p_cmZsparse) cholmod_l_free_sparse(&p_cmZsparse, cc);
if (p_Zsparse && *p_Zsparse) free(*p_Zsparse);
if (p_cmR) cholmod_l_free_sparse(&p_cmR, cc);
if (p_R && *p_R) free(*p_R);
if (p_cmH) cholmod_l_free_sparse(&p_cmH, cc);
if (p_H && *p_H) free(*p_H);
}
return rank;
}
#ifdef DEBUG_PRINT_MATRIX_VECTOR
void print_matrix_cmd(char* str, cholmod_dense* Md)
{
printf("%s :\n", str);
int i, j;
for (i = 0; i < Md->nrow; i++)
{
for (j = 0; j < Md->ncol; j++)
printf("%.16f ", ((double*)Md->x)[i + Md->nrow * j]);
printf("\n");
}
printf("\n");
}
void print_vector_cmd(char* str, cholmod_dense* Md)
{
printf("%s :\n", str);
int i;
for (i = 0; i < Md->nrow; i++)
{
printf("%.16f ", ((double*)Md->x)[i]);
}
printf("\n");
}
void print_matrix_cms(char* str, cholmod_sparse* Ms, cholmod_common* c)
{
cholmod_dense* Md = cholmod_l_sparse_to_dense(Ms, c);
print_matrix_cmd(str, Md);
cholmod_l_free_dense(&Md, c);
}
void print_matrix_s(char* str, splm_ccsm* Ms)
{
printf("%s :\n", str);
SuiteSparse_long i, j, idx;
for (i = 0; i < Ms->nr; i++)
{
for (j = 0; j < Ms->nc; j++)
{
idx = splm_ccsm_elmidx(Ms, i, j);
printf("%.16f ", idx >= 0 ? Ms->val[idx] : 0);
}
printf("\n");
}
printf("\n");
}
#endif
/* finite difference approximation to the Jacobian of func
* using either forward or central differences.
* The structure of the Jacobian is assumed known.
* Uses the strategy described in Nocedal-Wright, ch. 7, pp. 169.
*/
static SuiteSparse_long splm_intern_fdif_jac(
double* p, /* input: current parameter estimate, nvarsx1 */
splm_ccsm* jac, /* output: CCS array storing approximate Jacobian, nobsxnvars */
int nobs,
int mvars,
double epsfcn,
double mindeltax,
int (*func)(int n_obs, int m_vars, double* x, double* fvec, int* iflag),
int forw, /* 1: forward differences, 0: central differences */
int *nfeval, /* number of func evaluations for computing Jacobian */
int *iflag)
{
register SuiteSparse_long i, j, jj, k;
register double d;
SuiteSparse_long ii, m, * jcol = 0, * varlist = 0, * coldone = 0;
SuiteSparse_long* vidxs = 0, * ridxs = 0;
double* tmpd = 0;
double* p0 = 0, * hx = 0, * hxx = 0;
const double eps = sqrt(fmax(epsfcn, MACHEP));
double scl;
SuiteSparse_long ret = -1;
if (0 == (p0 = (double*)malloc(mvars * sizeof(double)))) goto freeall;
if (0 == (hx = (double*)malloc(nobs * sizeof(double)))) goto freeall;
if (0 == (hxx = (double*)malloc(nobs * sizeof(double)))) goto freeall;
if (0 == (jcol = (SuiteSparse_long*)malloc(nobs * sizeof(SuiteSparse_long)))) goto freeall; /* keeps track of measurements influenced by the set of variables currently in "varlist" below */
for (i = 0; i < nobs; ++i) jcol[i] = -1;
k = splm_ccsm_col_maxnelms(jac);
if (0 == (vidxs = (SuiteSparse_long*)malloc(2 * k * sizeof(SuiteSparse_long)))) goto freeall;
ridxs = vidxs + k;
if (0 == (varlist = (SuiteSparse_long*)malloc(mvars * sizeof(SuiteSparse_long)))) goto freeall; /* stores indices of J's columns which are computed with the same "func" call */
if (0 == (coldone = (SuiteSparse_long*)malloc(mvars * sizeof(SuiteSparse_long)))) goto freeall; /* keeps track of J's columns which have been already computed */
memset(coldone, 0, mvars * sizeof(SuiteSparse_long)); /* initialize to zero */
if (0 == (tmpd = (double*)malloc(mvars * sizeof(double)))) goto freeall;
if (forw)
{
(*func)(nobs, mvars, p, hx, iflag); ++(*nfeval); // hx=f(p)
}
for (j = 0; j < mvars; ++j)
{
p0[j] = p[j];
}
for (j = 0; j < mvars; ++j)
{
if (coldone[j]) continue; /* column j already computed */
k = splm_ccsm_col_elmidxs(jac, j, vidxs, ridxs);
for (i = 0; i < k; ++i) jcol[ridxs[i]] = j;
varlist[0] = j; m = 1; coldone[j] = 1;
for (jj = j + 1; jj < mvars; ++jj)
{
if (coldone[jj]) continue; /* column jj already computed */
k = splm_ccsm_col_elmidxs(jac, jj, vidxs, ridxs);
for (i = 0; i < k; ++i)
if (jcol[ridxs[i]] != -1) goto nextjj;
if (0 == k) { coldone[jj] = 1; continue; } /* all zeros column, ignore */
/* column jj does not clash with previously considered ones, mark it */
for (i = 0; i < k; ++i) jcol[ridxs[i]] = jj;
varlist[m++] = jj; coldone[jj] = 1;
nextjj:
continue;
}
for (k = 0; k < m; ++k)
{
/* determine d=max(SPLM_DELTA_SCALE*|p[varlist[k]]|, delta), see HZ */
d = eps * p[varlist[k]]; // force evaluation
d = fabs(d);
if (d < mindeltax) d = mindeltax;
if (!forw) d *= 0.5;
if (0 == d) d = eps;
tmpd[varlist[k]] = d;
p[varlist[k]] += d;
}
(*func)(nobs, mvars, p, hxx, iflag); ++(*nfeval); // hxx=f(p+d)
if (forw)
{
for (k = 0; k < m; ++k)
p[varlist[k]] = p0[varlist[k]]; /* restore */
scl = 1.0;
}
else
{ // central
for (k = 0; k < m; ++k)
p[varlist[k]] -= 2 * tmpd[varlist[k]];
(*func)(nobs, mvars, p, hx, iflag); ++(*nfeval); // hx=f(p-d)
for (k = 0; k < m; ++k)
p[varlist[k]] = p0[varlist[k]]; /* restore */
scl = 0.5; // 1./2.
}
for (k = 0; k < m; ++k)
{
d = tmpd[varlist[k]];
d = scl / d; /* invert so that divisions can be carried out faster as multiplications */
jj = splm_ccsm_col_elmidxs(jac, varlist[k], vidxs, ridxs);
for (i = 0; i < jj; ++i)
{
ii = ridxs[i];
jac->val[vidxs[i]] = (hxx[ii] - hx[ii]) * d;
jcol[ii] = -1; /* restore */
}
}
}
ret = 0;
freeall:
if (p0) free(p0);
if (hx) free(hx);
if (hxx) free(hxx);
if (tmpd) free(tmpd);
if (coldone) free(coldone);
if (varlist) free(varlist);
if (vidxs) free(vidxs);
if (jcol) free(jcol);
return ret;
}
cholmod_dense* dense_wrapper
(
cholmod_dense* X,
SuiteSparse_long nrow,
SuiteSparse_long ncol,
double* Xx
)
{
X->xtype = CHOLMOD_REAL;
X->nrow = nrow;
X->ncol = ncol;
X->d = nrow; // leading dimension = nrow
X->nzmax = nrow * ncol;
X->x = Xx;
X->z = 0; // ZOMPLEX case not supported
X->dtype = CHOLMOD_DOUBLE;
return X;
}
// =============================================================================
// === Rsolve ==================================================================
// =============================================================================
// Solve R*Y = X for Y and write Y back to X
SuiteSparse_long Rsolve
(
// R is at n-by-n or bigger, upper triangular with zero-free diagonal
// Use only the left, upper n-by-n submatrix of R
// X has memory for at least n doubles.
SuiteSparse_long n,
cholmod_sparse* R,
double* X, // X is n-by-nx, leading dimension n, overwritten with solution
SuiteSparse_long nx, // nx==1 if X is a column vector
cholmod_common* cc
)
{
SuiteSparse_long* Rp = (SuiteSparse_long*)R->p;
SuiteSparse_long* Ri = (SuiteSparse_long*)R->i;
double* Rx = (double*)R->x;
double rjj;
SuiteSparse_long k, j, p;
if (!R->packed || !R->sorted || R->dtype != CHOLMOD_DOUBLE || R->xtype != CHOLMOD_REAL)
{
return -1;
}
if ((SuiteSparse_long)R->ncol < n || (SuiteSparse_long)R->nrow < n)
{
return -1;
}
// check the diagonal
for (j = 0; j < n; j++)
{
if (Rp[j + 1] <= Rp[j] || Ri[Rp[j + 1] - 1] != j)
{
// Rsolve: R not upper triangular with zero-free diagonal
return -1;
}
}
// do the backsolve
for (k = 0; k < nx; k++)
{
for (j = n; j > n; )
{
X[--j] = ((double)0);
}
for (j = n; 0 != j; )
{
rjj = Rx[Rp[j] - 1];
if (((double)0) == rjj)
{
// Rsolve: R has an explicit zero on the diagonal
return -1;
}
X[--j] /= rjj;
for (p = Rp[j]; p + 1 < Rp[j + 1]; p++)
{
X[Ri[p]] -= Rx[p] * X[j];
}
}
X += n;
}
return 0;
}
// =============================================================================
// === RTsolve =================================================================
// =============================================================================
// Solve transpose(R)*Y = X for Y and write Y back to X
SuiteSparse_long RTsolve
(
// R is at n-by-n or bigger, upper triangular with zero-free diagonal
// Use only the left, upper n-by-n submatrix of R
// X has memory for at least n doubles.
SuiteSparse_long n,
cholmod_sparse* R,
double* X, // X is n-by-nx, leading dimension n, overwritten with solution
SuiteSparse_long nx, // nx==1 if X is a column vector
cholmod_common* cc
)
{
SuiteSparse_long i, j, k, l, p;
SuiteSparse_long* rowidx = R->i, * colptr = R->p;
SuiteSparse_long* colidx, * rowptr;
SuiteSparse_long* rowcounts; // counters for the number of nonzeros in each row
double* val;
double rjj;
if (!R->packed || !R->sorted || R->dtype != CHOLMOD_DOUBLE || R->xtype != CHOLMOD_REAL)
{
return -1;
}
if ((SuiteSparse_long)R->ncol < n || (SuiteSparse_long)R->nrow < n)
{
return -1;
}
// check the diagonal
for (j = 0; j < n; j++)
{
if (colptr[j + 1] <= colptr[j] || rowidx[colptr[j + 1] - 1] != j)
{
// RTsolve: R not upper triangular with zero-free diagonal
return -1;
}
}
if (0 == (rowcounts = (SuiteSparse_long*)malloc(n * sizeof(SuiteSparse_long))))
{
return -1;
}
if (0 == (val = (double*)malloc(R->nzmax * sizeof(double))))
{
free(rowcounts);
return -1;
}
if (0 == (colidx = (SuiteSparse_long*)malloc(R->nzmax * sizeof(SuiteSparse_long))))
{
free(rowcounts);
free(val);
return -1;
}
if (0 == (rowptr = (SuiteSparse_long*)malloc((n + 1) * sizeof(SuiteSparse_long))))
{
free(rowcounts);
free(val);
free(colidx);
return -1;
}
/* 1st pass: count #nonzeros in each row */
for (j = 0; j < n; j++)
rowcounts[j] = 0;
for (i = colptr[n]; i-- > 0; )
++(rowcounts[rowidx[i]]);
/* 2nd pass: copy every nonzero to its right position into the CRS structure */
for (i = k = 0; i < n; ++i) {
rowptr[i] = k;
k += rowcounts[i];
rowcounts[i] = 0; // clear to avoid separate loop below
}
rowptr[n] = R->nzmax;
/* rowcounts[i] will count the #nonzeros in row i seen before the current column */