1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 package org.apache.commons.math4.legacy.fitting.leastsquares;
18
19 import java.util.Arrays;
20
21 import org.apache.commons.math4.legacy.exception.ConvergenceException;
22 import org.apache.commons.math4.legacy.exception.util.LocalizedFormats;
23 import org.apache.commons.math4.legacy.fitting.leastsquares.LeastSquaresProblem.Evaluation;
24 import org.apache.commons.math4.legacy.linear.ArrayRealVector;
25 import org.apache.commons.math4.legacy.linear.RealMatrix;
26 import org.apache.commons.math4.legacy.optim.ConvergenceChecker;
27 import org.apache.commons.math4.core.jdkmath.JdkMath;
28 import org.apache.commons.math4.legacy.core.IntegerSequence;
29 import org.apache.commons.numbers.core.Precision;
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 public class LevenbergMarquardtOptimizer implements LeastSquaresOptimizer {
112
113
114 private static final double TWO_EPS = 2 * Precision.EPSILON;
115
116
117
118 private final double initialStepBoundFactor;
119
120 private final double costRelativeTolerance;
121
122 private final double parRelativeTolerance;
123
124
125 private final double orthoTolerance;
126
127 private final double qrRankingThreshold;
128
129
130
131
132
133
134
135
136
137
138
139
140 public LevenbergMarquardtOptimizer() {
141 this(100, 1e-10, 1e-10, 1e-10, Precision.SAFE_MIN);
142 }
143
144
145
146
147
148
149
150
151
152
153
154
155 public LevenbergMarquardtOptimizer(
156 final double initialStepBoundFactor,
157 final double costRelativeTolerance,
158 final double parRelativeTolerance,
159 final double orthoTolerance,
160 final double qrRankingThreshold) {
161 this.initialStepBoundFactor = initialStepBoundFactor;
162 this.costRelativeTolerance = costRelativeTolerance;
163 this.parRelativeTolerance = parRelativeTolerance;
164 this.orthoTolerance = orthoTolerance;
165 this.qrRankingThreshold = qrRankingThreshold;
166 }
167
168
169
170
171
172
173
174
175
176
177
178 public LevenbergMarquardtOptimizer withInitialStepBoundFactor(double newInitialStepBoundFactor) {
179 return new LevenbergMarquardtOptimizer(
180 newInitialStepBoundFactor,
181 costRelativeTolerance,
182 parRelativeTolerance,
183 orthoTolerance,
184 qrRankingThreshold);
185 }
186
187
188
189
190
191 public LevenbergMarquardtOptimizer withCostRelativeTolerance(double newCostRelativeTolerance) {
192 return new LevenbergMarquardtOptimizer(
193 initialStepBoundFactor,
194 newCostRelativeTolerance,
195 parRelativeTolerance,
196 orthoTolerance,
197 qrRankingThreshold);
198 }
199
200
201
202
203
204
205 public LevenbergMarquardtOptimizer withParameterRelativeTolerance(double newParRelativeTolerance) {
206 return new LevenbergMarquardtOptimizer(
207 initialStepBoundFactor,
208 costRelativeTolerance,
209 newParRelativeTolerance,
210 orthoTolerance,
211 qrRankingThreshold);
212 }
213
214
215
216
217
218
219
220
221 public LevenbergMarquardtOptimizer withOrthoTolerance(double newOrthoTolerance) {
222 return new LevenbergMarquardtOptimizer(
223 initialStepBoundFactor,
224 costRelativeTolerance,
225 parRelativeTolerance,
226 newOrthoTolerance,
227 qrRankingThreshold);
228 }
229
230
231
232
233
234
235
236
237 public LevenbergMarquardtOptimizer withRankingThreshold(double newQRRankingThreshold) {
238 return new LevenbergMarquardtOptimizer(
239 initialStepBoundFactor,
240 costRelativeTolerance,
241 parRelativeTolerance,
242 orthoTolerance,
243 newQRRankingThreshold);
244 }
245
246
247
248
249
250
251
252 public double getInitialStepBoundFactor() {
253 return initialStepBoundFactor;
254 }
255
256
257
258
259
260
261
262 public double getCostRelativeTolerance() {
263 return costRelativeTolerance;
264 }
265
266
267
268
269
270
271
272 public double getParameterRelativeTolerance() {
273 return parRelativeTolerance;
274 }
275
276
277
278
279
280
281
282 public double getOrthoTolerance() {
283 return orthoTolerance;
284 }
285
286
287
288
289
290
291
292 public double getRankingThreshold() {
293 return qrRankingThreshold;
294 }
295
296
297 @Override
298 public Optimum optimize(final LeastSquaresProblem problem) {
299
300 final int nR = problem.getObservationSize();
301 final int nC = problem.getParameterSize();
302
303 final IntegerSequence.Incrementor iterationCounter = problem.getIterationCounter();
304 final IntegerSequence.Incrementor evaluationCounter = problem.getEvaluationCounter();
305
306 final ConvergenceChecker<Evaluation> checker = problem.getConvergenceChecker();
307
308
309 final int solvedCols = JdkMath.min(nR, nC);
310
311 double[] lmDir = new double[nC];
312
313 double lmPar = 0;
314
315
316 double delta = 0;
317 double xNorm = 0;
318 double[] diag = new double[nC];
319 double[] oldX = new double[nC];
320 double[] oldRes = new double[nR];
321 double[] qtf = new double[nR];
322 double[] work1 = new double[nC];
323 double[] work2 = new double[nC];
324 double[] work3 = new double[nC];
325
326
327
328 evaluationCounter.increment();
329
330 Evaluation current = problem.evaluate(problem.getStart());
331 double[] currentResiduals = current.getResiduals().toArray();
332 double currentCost = current.getCost();
333 double[] currentPoint = current.getPoint().toArray();
334
335
336 boolean firstIteration = true;
337 while (true) {
338 iterationCounter.increment();
339
340 final Evaluation previous = current;
341
342
343 final InternalData internalData
344 = qrDecomposition(current.getJacobian(), solvedCols);
345 final double[][] weightedJacobian = internalData.weightedJacobian;
346 final int[] permutation = internalData.permutation;
347 final double[] diagR = internalData.diagR;
348 final double[] jacNorm = internalData.jacNorm;
349
350
351 double[] weightedResidual = currentResiduals;
352 System.arraycopy(weightedResidual, 0, qtf, 0, nR);
353
354
355 qTy(qtf, internalData);
356
357
358
359 for (int k = 0; k < solvedCols; ++k) {
360 int pk = permutation[k];
361 weightedJacobian[k][pk] = diagR[pk];
362 }
363
364 if (firstIteration) {
365
366
367 xNorm = 0;
368 for (int k = 0; k < nC; ++k) {
369 double dk = jacNorm[k];
370 if (dk == 0) {
371 dk = 1.0;
372 }
373 double xk = dk * currentPoint[k];
374 xNorm += xk * xk;
375 diag[k] = dk;
376 }
377 xNorm = JdkMath.sqrt(xNorm);
378
379
380 delta = (xNorm == 0) ? initialStepBoundFactor : (initialStepBoundFactor * xNorm);
381 }
382
383
384 double maxCosine = 0;
385 if (currentCost != 0) {
386 for (int j = 0; j < solvedCols; ++j) {
387 int pj = permutation[j];
388 double s = jacNorm[pj];
389 if (s != 0) {
390 double sum = 0;
391 for (int i = 0; i <= j; ++i) {
392 sum += weightedJacobian[i][pj] * qtf[i];
393 }
394 maxCosine = JdkMath.max(maxCosine, JdkMath.abs(sum) / (s * currentCost));
395 }
396 }
397 }
398 if (maxCosine <= orthoTolerance) {
399
400 return new OptimumImpl(
401 current,
402 evaluationCounter.getCount(),
403 iterationCounter.getCount());
404 }
405
406
407 for (int j = 0; j < nC; ++j) {
408 diag[j] = JdkMath.max(diag[j], jacNorm[j]);
409 }
410
411
412 for (double ratio = 0; ratio < 1.0e-4;) {
413
414
415 for (int j = 0; j < solvedCols; ++j) {
416 int pj = permutation[j];
417 oldX[pj] = currentPoint[pj];
418 }
419 final double previousCost = currentCost;
420 double[] tmpVec = weightedResidual;
421 weightedResidual = oldRes;
422 oldRes = tmpVec;
423
424
425 lmPar = determineLMParameter(qtf, delta, diag,
426 internalData, solvedCols,
427 work1, work2, work3, lmDir, lmPar);
428
429
430 double lmNorm = 0;
431 for (int j = 0; j < solvedCols; ++j) {
432 int pj = permutation[j];
433 lmDir[pj] = -lmDir[pj];
434 currentPoint[pj] = oldX[pj] + lmDir[pj];
435 double s = diag[pj] * lmDir[pj];
436 lmNorm += s * s;
437 }
438 lmNorm = JdkMath.sqrt(lmNorm);
439
440 if (firstIteration) {
441 delta = JdkMath.min(delta, lmNorm);
442 }
443
444
445 evaluationCounter.increment();
446 current = problem.evaluate(new ArrayRealVector(currentPoint));
447 currentResiduals = current.getResiduals().toArray();
448 currentCost = current.getCost();
449 currentPoint = current.getPoint().toArray();
450
451
452 double actRed = -1.0;
453 if (0.1 * currentCost < previousCost) {
454 double r = currentCost / previousCost;
455 actRed = 1.0 - r * r;
456 }
457
458
459
460 for (int j = 0; j < solvedCols; ++j) {
461 int pj = permutation[j];
462 double dirJ = lmDir[pj];
463 work1[j] = 0;
464 for (int i = 0; i <= j; ++i) {
465 work1[i] += weightedJacobian[i][pj] * dirJ;
466 }
467 }
468 double coeff1 = 0;
469 for (int j = 0; j < solvedCols; ++j) {
470 coeff1 += work1[j] * work1[j];
471 }
472 double pc2 = previousCost * previousCost;
473 coeff1 /= pc2;
474 double coeff2 = lmPar * lmNorm * lmNorm / pc2;
475 double preRed = coeff1 + 2 * coeff2;
476 double dirDer = -(coeff1 + coeff2);
477
478
479 ratio = (preRed == 0) ? 0 : (actRed / preRed);
480
481
482 if (ratio <= 0.25) {
483 double tmp =
484 (actRed < 0) ? (0.5 * dirDer / (dirDer + 0.5 * actRed)) : 0.5;
485 if (0.1 * currentCost >= previousCost || tmp < 0.1) {
486 tmp = 0.1;
487 }
488 delta = tmp * JdkMath.min(delta, 10.0 * lmNorm);
489 lmPar /= tmp;
490 } else if (lmPar == 0 || ratio >= 0.75) {
491 delta = 2 * lmNorm;
492 lmPar *= 0.5;
493 }
494
495
496 if (ratio >= 1.0e-4) {
497
498 firstIteration = false;
499 xNorm = 0;
500 for (int k = 0; k < nC; ++k) {
501 double xK = diag[k] * currentPoint[k];
502 xNorm += xK * xK;
503 }
504 xNorm = JdkMath.sqrt(xNorm);
505
506
507 if (checker != null && checker.converged(iterationCounter.getCount(), previous, current)) {
508 return new OptimumImpl(current, evaluationCounter.getCount(), iterationCounter.getCount());
509 }
510 } else {
511
512 currentCost = previousCost;
513 for (int j = 0; j < solvedCols; ++j) {
514 int pj = permutation[j];
515 currentPoint[pj] = oldX[pj];
516 }
517 tmpVec = weightedResidual;
518 weightedResidual = oldRes;
519 oldRes = tmpVec;
520
521 current = previous;
522 }
523
524
525 if ((JdkMath.abs(actRed) <= costRelativeTolerance &&
526 preRed <= costRelativeTolerance &&
527 ratio <= 2.0) ||
528 delta <= parRelativeTolerance * xNorm) {
529 return new OptimumImpl(current, evaluationCounter.getCount(), iterationCounter.getCount());
530 }
531
532
533 if (JdkMath.abs(actRed) <= TWO_EPS &&
534 preRed <= TWO_EPS &&
535 ratio <= 2.0) {
536 throw new ConvergenceException(LocalizedFormats.TOO_SMALL_COST_RELATIVE_TOLERANCE,
537 costRelativeTolerance);
538 } else if (delta <= TWO_EPS * xNorm) {
539 throw new ConvergenceException(LocalizedFormats.TOO_SMALL_PARAMETERS_RELATIVE_TOLERANCE,
540 parRelativeTolerance);
541 } else if (maxCosine <= TWO_EPS) {
542 throw new ConvergenceException(LocalizedFormats.TOO_SMALL_ORTHOGONALITY_TOLERANCE,
543 orthoTolerance);
544 }
545 }
546 }
547 }
548
549
550
551
552
553
554
555 private static final class InternalData {
556
557 private final double[][] weightedJacobian;
558
559 private final int[] permutation;
560
561 private final int rank;
562
563 private final double[] diagR;
564
565 private final double[] jacNorm;
566
567 private final double[] beta;
568
569
570
571
572
573
574
575
576
577 InternalData(double[][] weightedJacobian,
578 int[] permutation,
579 int rank,
580 double[] diagR,
581 double[] jacNorm,
582 double[] beta) {
583 this.weightedJacobian = weightedJacobian;
584 this.permutation = permutation;
585 this.rank = rank;
586 this.diagR = diagR;
587 this.jacNorm = jacNorm;
588 this.beta = beta;
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 private double determineLMParameter(double[] qy, double delta, double[] diag,
621 InternalData internalData, int solvedCols,
622 double[] work1, double[] work2, double[] work3,
623 double[] lmDir, double lmPar) {
624 final double[][] weightedJacobian = internalData.weightedJacobian;
625 final int[] permutation = internalData.permutation;
626 final int rank = internalData.rank;
627 final double[] diagR = internalData.diagR;
628
629 final int nC = weightedJacobian[0].length;
630
631
632
633 for (int j = 0; j < rank; ++j) {
634 lmDir[permutation[j]] = qy[j];
635 }
636 for (int j = rank; j < nC; ++j) {
637 lmDir[permutation[j]] = 0;
638 }
639 for (int k = rank - 1; k >= 0; --k) {
640 int pk = permutation[k];
641 double ypk = lmDir[pk] / diagR[pk];
642 for (int i = 0; i < k; ++i) {
643 lmDir[permutation[i]] -= ypk * weightedJacobian[i][pk];
644 }
645 lmDir[pk] = ypk;
646 }
647
648
649
650 double dxNorm = 0;
651 for (int j = 0; j < solvedCols; ++j) {
652 int pj = permutation[j];
653 double s = diag[pj] * lmDir[pj];
654 work1[pj] = s;
655 dxNorm += s * s;
656 }
657 dxNorm = JdkMath.sqrt(dxNorm);
658 double fp = dxNorm - delta;
659 if (fp <= 0.1 * delta) {
660 lmPar = 0;
661 return lmPar;
662 }
663
664
665
666
667 double sum2;
668 double parl = 0;
669 if (rank == solvedCols) {
670 for (int j = 0; j < solvedCols; ++j) {
671 int pj = permutation[j];
672 work1[pj] *= diag[pj] / dxNorm;
673 }
674 sum2 = 0;
675 for (int j = 0; j < solvedCols; ++j) {
676 int pj = permutation[j];
677 double sum = 0;
678 for (int i = 0; i < j; ++i) {
679 sum += weightedJacobian[i][pj] * work1[permutation[i]];
680 }
681 double s = (work1[pj] - sum) / diagR[pj];
682 work1[pj] = s;
683 sum2 += s * s;
684 }
685 parl = fp / (delta * sum2);
686 }
687
688
689 sum2 = 0;
690 for (int j = 0; j < solvedCols; ++j) {
691 int pj = permutation[j];
692 double sum = 0;
693 for (int i = 0; i <= j; ++i) {
694 sum += weightedJacobian[i][pj] * qy[i];
695 }
696 sum /= diag[pj];
697 sum2 += sum * sum;
698 }
699 double gNorm = JdkMath.sqrt(sum2);
700 double paru = gNorm / delta;
701 if (paru == 0) {
702 paru = Precision.SAFE_MIN / JdkMath.min(delta, 0.1);
703 }
704
705
706
707 lmPar = JdkMath.min(paru, JdkMath.max(lmPar, parl));
708 if (lmPar == 0) {
709 lmPar = gNorm / dxNorm;
710 }
711
712 for (int countdown = 10; countdown >= 0; --countdown) {
713
714
715 if (lmPar == 0) {
716 lmPar = JdkMath.max(Precision.SAFE_MIN, 0.001 * paru);
717 }
718 double sPar = JdkMath.sqrt(lmPar);
719 for (int j = 0; j < solvedCols; ++j) {
720 int pj = permutation[j];
721 work1[pj] = sPar * diag[pj];
722 }
723 determineLMDirection(qy, work1, work2, internalData, solvedCols, work3, lmDir);
724
725 dxNorm = 0;
726 for (int j = 0; j < solvedCols; ++j) {
727 int pj = permutation[j];
728 double s = diag[pj] * lmDir[pj];
729 work3[pj] = s;
730 dxNorm += s * s;
731 }
732 dxNorm = JdkMath.sqrt(dxNorm);
733 double previousFP = fp;
734 fp = dxNorm - delta;
735
736
737
738 if (JdkMath.abs(fp) <= 0.1 * delta ||
739 (parl == 0 &&
740 fp <= previousFP &&
741 previousFP < 0)) {
742 return lmPar;
743 }
744
745
746 for (int j = 0; j < solvedCols; ++j) {
747 int pj = permutation[j];
748 work1[pj] = work3[pj] * diag[pj] / dxNorm;
749 }
750 for (int j = 0; j < solvedCols; ++j) {
751 int pj = permutation[j];
752 work1[pj] /= work2[j];
753 double tmp = work1[pj];
754 for (int i = j + 1; i < solvedCols; ++i) {
755 work1[permutation[i]] -= weightedJacobian[i][pj] * tmp;
756 }
757 }
758 sum2 = 0;
759 for (int j = 0; j < solvedCols; ++j) {
760 double s = work1[permutation[j]];
761 sum2 += s * s;
762 }
763 double correction = fp / (delta * sum2);
764
765
766 if (fp > 0) {
767 parl = JdkMath.max(parl, lmPar);
768 } else if (fp < 0) {
769 paru = JdkMath.min(paru, lmPar);
770 }
771
772
773 lmPar = JdkMath.max(parl, lmPar + correction);
774 }
775
776 return lmPar;
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 private void determineLMDirection(double[] qy, double[] diag,
803 double[] lmDiag,
804 InternalData internalData,
805 int solvedCols,
806 double[] work,
807 double[] lmDir) {
808 final int[] permutation = internalData.permutation;
809 final double[][] weightedJacobian = internalData.weightedJacobian;
810 final double[] diagR = internalData.diagR;
811
812
813
814 for (int j = 0; j < solvedCols; ++j) {
815 int pj = permutation[j];
816 for (int i = j + 1; i < solvedCols; ++i) {
817 weightedJacobian[i][pj] = weightedJacobian[j][permutation[i]];
818 }
819 lmDir[j] = diagR[pj];
820 work[j] = qy[j];
821 }
822
823
824 for (int j = 0; j < solvedCols; ++j) {
825
826
827
828 int pj = permutation[j];
829 double dpj = diag[pj];
830 if (dpj != 0) {
831 Arrays.fill(lmDiag, j + 1, lmDiag.length, 0);
832 }
833 lmDiag[j] = dpj;
834
835
836
837
838 double qtbpj = 0;
839 for (int k = j; k < solvedCols; ++k) {
840 int pk = permutation[k];
841
842
843
844 if (lmDiag[k] != 0) {
845
846 final double sin;
847 final double cos;
848 double rkk = weightedJacobian[k][pk];
849 if (JdkMath.abs(rkk) < JdkMath.abs(lmDiag[k])) {
850 final double cotan = rkk / lmDiag[k];
851 sin = 1.0 / JdkMath.sqrt(1.0 + cotan * cotan);
852 cos = sin * cotan;
853 } else {
854 final double tan = lmDiag[k] / rkk;
855 cos = 1.0 / JdkMath.sqrt(1.0 + tan * tan);
856 sin = cos * tan;
857 }
858
859
860
861 weightedJacobian[k][pk] = cos * rkk + sin * lmDiag[k];
862 final double temp = cos * work[k] + sin * qtbpj;
863 qtbpj = -sin * work[k] + cos * qtbpj;
864 work[k] = temp;
865
866
867 for (int i = k + 1; i < solvedCols; ++i) {
868 double rik = weightedJacobian[i][pk];
869 final double temp2 = cos * rik + sin * lmDiag[i];
870 lmDiag[i] = -sin * rik + cos * lmDiag[i];
871 weightedJacobian[i][pk] = temp2;
872 }
873 }
874 }
875
876
877
878 lmDiag[j] = weightedJacobian[j][permutation[j]];
879 weightedJacobian[j][permutation[j]] = lmDir[j];
880 }
881
882
883
884 int nSing = solvedCols;
885 for (int j = 0; j < solvedCols; ++j) {
886 if (lmDiag[j] == 0 && nSing == solvedCols) {
887 nSing = j;
888 }
889 if (nSing < solvedCols) {
890 work[j] = 0;
891 }
892 }
893 if (nSing > 0) {
894 for (int j = nSing - 1; j >= 0; --j) {
895 int pj = permutation[j];
896 double sum = 0;
897 for (int i = j + 1; i < nSing; ++i) {
898 sum += weightedJacobian[i][pj] * work[i];
899 }
900 work[j] = (work[j] - sum) / lmDiag[j];
901 }
902 }
903
904
905 for (int j = 0; j < lmDir.length; ++j) {
906 lmDir[permutation[j]] = work[j];
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 private InternalData qrDecomposition(RealMatrix jacobian,
937 int solvedCols) throws ConvergenceException {
938
939
940 final double[][] weightedJacobian = jacobian.scalarMultiply(-1).getData();
941
942 final int nR = weightedJacobian.length;
943 final int nC = weightedJacobian[0].length;
944
945 final int[] permutation = new int[nC];
946 final double[] diagR = new double[nC];
947 final double[] jacNorm = new double[nC];
948 final double[] beta = new double[nC];
949
950
951 for (int k = 0; k < nC; ++k) {
952 permutation[k] = k;
953 double norm2 = 0;
954 for (int i = 0; i < nR; ++i) {
955 double akk = weightedJacobian[i][k];
956 norm2 += akk * akk;
957 }
958 jacNorm[k] = JdkMath.sqrt(norm2);
959 }
960
961
962 for (int k = 0; k < nC; ++k) {
963
964
965 int nextColumn = -1;
966 double ak2 = Double.NEGATIVE_INFINITY;
967 for (int i = k; i < nC; ++i) {
968 double norm2 = 0;
969 for (int j = k; j < nR; ++j) {
970 double aki = weightedJacobian[j][permutation[i]];
971 norm2 += aki * aki;
972 }
973 if (Double.isInfinite(norm2) || Double.isNaN(norm2)) {
974 throw new ConvergenceException(LocalizedFormats.UNABLE_TO_PERFORM_QR_DECOMPOSITION_ON_JACOBIAN,
975 nR, nC);
976 }
977 if (norm2 > ak2) {
978 nextColumn = i;
979 ak2 = norm2;
980 }
981 }
982 if (ak2 <= qrRankingThreshold) {
983 return new InternalData(weightedJacobian, permutation, k, diagR, jacNorm, beta);
984 }
985 int pk = permutation[nextColumn];
986 permutation[nextColumn] = permutation[k];
987 permutation[k] = pk;
988
989
990 double akk = weightedJacobian[k][pk];
991 double alpha = (akk > 0) ? -JdkMath.sqrt(ak2) : JdkMath.sqrt(ak2);
992 double betak = 1.0 / (ak2 - akk * alpha);
993 beta[pk] = betak;
994
995
996 diagR[pk] = alpha;
997 weightedJacobian[k][pk] -= alpha;
998
999
1000 for (int dk = nC - 1 - k; dk > 0; --dk) {
1001 double gamma = 0;
1002 for (int j = k; j < nR; ++j) {
1003 gamma += weightedJacobian[j][pk] * weightedJacobian[j][permutation[k + dk]];
1004 }
1005 gamma *= betak;
1006 for (int j = k; j < nR; ++j) {
1007 weightedJacobian[j][permutation[k + dk]] -= gamma * weightedJacobian[j][pk];
1008 }
1009 }
1010 }
1011
1012 return new InternalData(weightedJacobian, permutation, solvedCols, diagR, jacNorm, beta);
1013 }
1014
1015
1016
1017
1018
1019
1020
1021 private void qTy(double[] y,
1022 InternalData internalData) {
1023 final double[][] weightedJacobian = internalData.weightedJacobian;
1024 final int[] permutation = internalData.permutation;
1025 final double[] beta = internalData.beta;
1026
1027 final int nR = weightedJacobian.length;
1028 final int nC = weightedJacobian[0].length;
1029
1030 for (int k = 0; k < nC; ++k) {
1031 int pk = permutation[k];
1032 double gamma = 0;
1033 for (int i = k; i < nR; ++i) {
1034 gamma += weightedJacobian[i][pk] * y[i];
1035 }
1036 gamma *= beta[pk];
1037 for (int i = k; i < nR; ++i) {
1038 y[i] -= gamma * weightedJacobian[i][pk];
1039 }
1040 }
1041 }
1042 }