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