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 initialize();
188
189 assert(f && (int)f->size == ncon);
190 assert(r && (int)r->size == ncon);
191 assert(Fetaxi && (int)Fetaxi->size1 == ncon && (int)Fetaxi->size2 == npar);
192 assert(S && (int)S->size1 == ncon && (int)S->size2 == ncon);
193 assert(Sinv && (int)Sinv->size1 == ncon && (int)Sinv->size2 == ncon);
194 assert(nunm == 0 || (SinvFxi && (int)SinvFxi->size1 == ncon && (int)SinvFxi->size2 == nunm));
195 assert(SinvFeta && (int)SinvFeta->size1 == ncon && (int)SinvFeta->size2 == nmea);
196 assert(nunm == 0 || (W1 && (int)W1->size1 == nunm && (int)W1->size2 == nunm));
197 assert(G && (int)G->size1 == nmea && (int)G->size2 == nmea);
198 assert(nunm == 0 || (H && (int)H->size1 == nmea && (int)H->size2 == nunm));
199 assert(nunm == 0 || (HU && (int)HU->size1 == nmea && (int)HU->size2 == nunm));
200 assert(IGV && (int)IGV->size1 == nmea && (int)IGV->size2 == nmea);
201 assert(V && (int)V->size1 == npar && (int)V->size2 == npar);
202 assert(VLU && (int)VLU->size1 == nmea && (int)VLU->size2 == nmea);
203 assert(Vinv && (int)Vinv->size1 == nmea && (int)Vinv->size2 == nmea);
204 assert(Vnew && (int)Vnew->size1 == npar && (int)Vnew->size2 == npar);
205 assert(Minv && (int)Minv->size1 == npar && (int)Minv->size2 == npar);
206 assert(dxdt && (int)dxdt->size1 == npar && (int)dxdt->size2 == nmea);
207 assert(Vdxdt && (int)Vdxdt->size1 == nmea && (int)Vdxdt->size2 == npar);
208 assert(nunm == 0 || (dxi && (int)dxi->size == nunm));
209 assert(nunm == 0 || (Fxidxi && (int)Fxidxi->size == ncon));
210 assert(lambda && (int)lambda->size == ncon);
211 assert(FetaTlambda && (int)FetaTlambda->size == nmea);
212 assert(etaxi && (int)etaxi->size == npar);
213 assert(etasv && (int)etasv->size == npar);
214 assert(y && (int)y->size == nmea);
215 assert(y_eta && (int)y_eta->size == nmea);
216 assert(Vinvy_eta && (int)Vinvy_eta->size == nmea);
217 assert(FetaV && (int)FetaV->size1 == ncon && (int)FetaV->size2 == nmea);
219 assert(nunm == 0 || (
permU && (
int)
permU->size == nunm));
221
222
223 gsl_vector_view eta = gsl_vector_subvector(etaxi, 0, nmea);
224
225
226
227 B2DEBUG(11, "==== " << ncon << " " << nmea);
228
229 gsl_matrix_view Feta = gsl_matrix_submatrix(Fetaxi, 0, 0, ncon, nmea);
230
231 gsl_matrix_view Vetaeta = gsl_matrix_submatrix(V, 0, 0, nmea, nmea);
232
233 for (unsigned int i = 0; i < fitobjects.size(); ++i) {
234 for (int ilocal = 0; ilocal < fitobjects[i]->getNPar(); ++ilocal) {
235 if (!fitobjects[i]->isParamFixed(ilocal)) {
236 int iglobal = fitobjects[i]->getGlobalParNum(ilocal);
237 assert(iglobal >= 0 && iglobal < npar);
238 gsl_vector_set(etaxi, iglobal, fitobjects[i]->getParam(ilocal));
239 if (fitobjects[i]->isParamMeasured(ilocal)) {
240 assert(iglobal < nmea);
241 gsl_vector_set(y, iglobal, fitobjects[i]->getMParam(ilocal));
242 }
243 B2DEBUG(10, "etaxi[" << iglobal << "] = " << gsl_vector_get(etaxi, iglobal)
244 << " for jet " << i << " and ilocal = " << ilocal);
245 }
246 }
247 }
248
250 gsl_matrix_set_zero(Fetaxi);
251 for (int k = 0; k < ncon; k++) {
252 constraints[k]->getDerivatives(Fetaxi->size2, Fetaxi->block->data + k * Fetaxi->tda);
253 if (debug > 1) for (int j = 0; j < npar; j++)
254 if (gsl_matrix_get(Fetaxi, k, j) != 0)
255 B2INFO("1: Fetaxi[" << k << "][" << j << "] = " << gsl_matrix_get(Fetaxi, k, j));
256 }
257
258
259 double chinew = 0, chit = 0, chik = 0;
260 double alph = 1.;
261 nit = 0;
262
263 int nitmax = 200;
264 double chik0 = 100.;
265 double chit0 = 100.;
266 double dchikc = 1.0E-3;
267 double dchitc = 1.0E-4;
268 double dchikt = 1.0E-2;
269 double dchik = 1.05;
270 double chimxw = 10000.;
271 double almin = 0.05;
272
273
274 bool repeat = true;
275 bool scut = false;
276 bool calcerr = true;
277
278#ifndef FIT_TRACEOFF
280#endif
281
282
283
284 while (repeat) {
285 bool updatesuccess = true;
286
287
288 if (scut) {
289 gsl_vector_memcpy(etaxi, etasv);
290 updatesuccess = updateFitObjects(etaxi->block->data);
291 if (!updatesuccess) {
292 B2ERROR("OPALFitterGSL::fit: old parameters are garbage!");
293 return -1;
294 }
295
296
297 gsl_matrix_set_zero(Fetaxi);
298 for (int k = 0; k < ncon; ++k) {
299 constraints[k]->getDerivatives(Fetaxi->size2, Fetaxi->block->data + k * Fetaxi->tda);
300 }
301 if (debug > 1) debug_print(Fetaxi, "1: Fetaxi");
302 } else {
303 gsl_vector_memcpy(etasv, etaxi);
304 chik0 = chik;
305 chit0 = chit;
306 }
307
308
309 gsl_matrix_set_zero(V);
310
311 for (auto& fitobject : fitobjects) {
312 fitobject->addToGlobCov(V->block->data, V->tda);
313 }
314 if (debug > 1) debug_print(V, "V");
315
316 gsl_matrix_memcpy(VLU, &Vetaeta.matrix);
317
318
319
320 int signum;
321 int result;
322 result = gsl_linalg_LU_decomp(VLU,
permV, &signum);
323 B2DEBUG(11, "gsl_linalg_LU_decomp result=" << result);
324 if (debug > 3) debug_print(VLU, "VLU");
325
326 result = gsl_linalg_LU_invert(VLU,
permV, Vinv);
327 B2DEBUG(11, "gsl_linalg_LU_invert result=" << result);
328
329 if (debug > 2) debug_print(Vinv, "Vinv");
330
331
332 for (int k = 0; k < ncon; ++k) {
333 gsl_vector_set(f, k, constraints[k]->getValue());
334 }
335 if (debug > 1) debug_print(f, "f");
336
337
338 gsl_vector_memcpy(y_eta, y);
339 gsl_vector_sub(y_eta, &eta.vector);
340
341 gsl_vector_memcpy(r, f);
342
343 gsl_blas_dgemv(CblasNoTrans, 1, &Feta.matrix, y_eta, 1, r);
344
345 if (debug > 1) debug_print(r, "r");
346
347
348
349
350 B2DEBUG(12, "Creating FetaV");
351 gsl_blas_dsymm(CblasRight, CblasUpper, 1, &Vetaeta.matrix, &Feta.matrix, 0, FetaV);
352
353 B2DEBUG(12, "Creating S");
354 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, FetaV, &Feta.matrix, 0, S);
355
356 if (nunm > 0) {
357
358
359
360
361
362 gsl_matrix_view Fxi = gsl_matrix_submatrix(Fetaxi, 0, nmea, ncon, nunm);
363
364
365 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, &Fxi.matrix, &Fxi.matrix, 1, S);
366 }
367
368 if (debug > 1) debug_print(S, "S");
369
370
371
372
373 gsl_linalg_LU_decomp(S,
permS, &signum);
374 int inverr = gsl_linalg_LU_invert(S,
permS, Sinv);
375
376 if (inverr != 0) {
377 B2ERROR("S: gsl_linalg_LU_invert error " << inverr);
378 ierr = 7;
379 calcerr = false;
380 break;
381 }
382
383
384
385
386 gsl_blas_dsymv(CblasUpper, 1, Sinv, r, 0, lambda);
387
388
389
390 if (nunm > 0) {
391
392 gsl_matrix_view Fxi = gsl_matrix_submatrix(Fetaxi, 0, nmea, ncon, nunm);
393
394
395 gsl_blas_dsymm(CblasLeft, CblasUpper, 1, Sinv, &Fxi.matrix, 0, SinvFxi);
396
397 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1, &Fxi.matrix, SinvFxi, 0, W1);
398
399 if (debug > 1) {
400 debug_print(W1, "W1");
401
402 for (int i = 0; i < nunm; ++i) {
403 for (int j = 0; j < nunm; ++j) {
404 if (abs(gsl_matrix_get(W1, i, j) - gsl_matrix_get(W1, j, i)) > 1E-3 * abs(gsl_matrix_get(W1, i, j) + gsl_matrix_get(W1, j, i)))
405 B2INFO("W1[" << i << "][" << j << "] = " << gsl_matrix_get(W1, i, j)
406 << " W1[" << j << "][" << i << "] = " << gsl_matrix_get(W1, j, i)
407 << " => diff=" << abs(gsl_matrix_get(W1, i, j) - gsl_matrix_get(W1, j, i))
408 << " => tol=" << 1E-3 * abs(gsl_matrix_get(W1, i, j) + gsl_matrix_get(W1, j, i)));
409 }
410 }
411 }
412
413
414
415
416
417
418
419
420
421 B2DEBUG(11, "alph = " << alph);
422 if (debug > 1) {
423 debug_print(lambda, "lambda");
424 debug_print(&(Fxi.matrix), "Fxi");
425 }
426
427 gsl_blas_dgemv(CblasTrans, -alph, &Fxi.matrix, lambda, 0, dxi);
428
429 if (debug > 1) {
430 debug_print(dxi, "dxi0");
431 debug_print(W1, "W1");
432 }
433
434
435
436
437 gsl_linalg_cholesky_decomp(W1);
438 inverr = gsl_linalg_cholesky_svx(W1, dxi);
439
440 if (debug > 1) debug_print(dxi, "dxi1");
441
442
443 if (inverr != 0) {
444 B2ERROR("W1: gsl_linalg_cholesky_svx error " << inverr);
445 ierr = 8;
446 calcerr = false;
447 break;
448 }
449
450 if (debug > 1) debug_print(dxi, "dxi");
451
452
453
454 gsl_vector_view xi = gsl_vector_subvector(etaxi, nmea, nunm);
455
456
457 gsl_vector_add(&xi.vector, dxi);
458
459 }
460
461
462
463
464
465
466 if (nunm > 0) {
467
468 gsl_matrix_view Fxi = gsl_matrix_submatrix(Fetaxi, 0, nmea, ncon, nunm);
469
470 gsl_blas_dgemv(CblasNoTrans, 1, &Fxi.matrix, dxi, 0, Fxidxi);
471
472 gsl_blas_dsymv(CblasUpper, 1, Sinv, Fxidxi, 1, lambda);
473
474 }
475
476 if (debug > 1) debug_print(lambda, "lambda");
477
478
479
480 gsl_vector_memcpy(&eta.vector, y);
481
482 gsl_blas_dgemv(CblasTrans, 1, &Feta.matrix, lambda, 0, FetaTlambda);
483
484 gsl_blas_dsymv(CblasUpper, -1, &Vetaeta.matrix, FetaTlambda, 1, &eta.vector);
485
486
487 if (debug > 1) debug_print(&eta.vector, "updated eta");
488
489
490
491
492
493 updatesuccess = updateFitObjects(etaxi->block->data);
494
495 B2DEBUG(10, "After adjustment of all parameters:\n");
496 for (int k = 0; k < ncon; ++k) {
497 B2DEBUG(10, "Value of constraint " << k << " = " << constraints[k]->getValue());
498 }
499 gsl_matrix_set_zero(Fetaxi);
500 for (int k = 0; k < ncon; k++) {
501 constraints[k]->getDerivatives(Fetaxi->size2, Fetaxi->block->data + k * Fetaxi->tda);
502 }
503 if (debug > 1) debug_print(Fetaxi, "2: Fetaxi");
504
505
506
507
508
509 gsl_vector_memcpy(y_eta, y);
510 gsl_vector_sub(y_eta, &eta.vector);
511
512
513 gsl_blas_dsymv(CblasUpper, 1, Vinv, y_eta, 0, Vinvy_eta);
514
515 gsl_blas_ddot(y_eta, Vinvy_eta, &chit);
516
517 for (int i = 0; i < nmea; ++i)
518 for (int j = 0; j < nmea; ++j) {
519 double dchit = (gsl_vector_get(y_eta, i)) *
520 gsl_matrix_get(Vinv, i, j) *
521 (gsl_vector_get(y_eta, j));
522 if (dchit != 0)
523 B2DEBUG(11, "chit for i,j = " << i << " , " << j << " = "
524 << dchit);
525 }
526
527 chik = 0;
528 for (int k = 0; k < ncon; ++k) chik += std::abs(2 * gsl_vector_get(lambda, k) * constraints[k]->getValue());
529 chinew = chit + chik;
530
531
532 nit++;
533
534 bool sconv = (std::abs(chik - chik0) < dchikc)
535 && (std::abs(chit - chit0) < dchitc * chit)
536 && (chik < dchikt * chit);
537
538
539
540
541
542
543
544 double eps = 1E-6;
545 bool sconv2 = true;
546 for (int k = 0; sconv2 && (k < ncon); ++k)
547 sconv2 &= (std::abs(gsl_vector_get(f, k)) < eps);
548 if (sconv2)
549 B2DEBUG(10, "All constraints fulfilled to better than " << eps);
550
551 for (int j = 0; sconv2 && (j < npar); ++j)
552 sconv2 &= (std::abs(gsl_vector_get(etaxi, j) - gsl_vector_get(etasv, j)) < eps);
553 if (sconv2)
554 B2DEBUG(10, "All parameters stable to better than " << eps);
555 sconv |= sconv2;
556
557 bool sbad = (chik > dchik * chik0)
558 && (chik > dchikt * chit)
559 && (chik > chik0 + 1.E-10);
560
561 scut = false;
562
563 if (nit > nitmax) {
564
565 repeat = false;
566 calcerr = false;
567 ierr = 1;
568 } else if (sconv && updatesuccess) {
569
570 repeat = false;
571 calcerr = true;
572 ierr = 0;
573 } else if (nit > 2 && chinew > chimxw && updatesuccess) {
574
575 repeat = false;
576 calcerr = false;
577 ierr = 2;
578 } else if ((sbad && nit > 1) || !updatesuccess) {
579
580 if (alph == almin) {
581 repeat = true;
582 ierr = 3;
583 } else {
584 alph = std::max(almin, 0.5 * alph);
585 scut = true;
586 repeat = true;
587 ierr = 4;
588 }
589 } else {
590
591 alph = std::min(alph + 0.1, 1.);
592 repeat = true;
593 ierr = 5;
594 }
595
596 B2DEBUG(10, "======== NIT = " << nit << ", CHI2 = " << chinew
597 << ", ierr = " << ierr << ", alph=" << alph);
598
599 for (unsigned int i = 0; i < fitobjects.size(); ++i)
600 B2DEBUG(10, "fitobject " << i << ": " << *fitobjects[i]);
601
602#ifndef FIT_TRACEOFF
603 if (tracer) tracer->
step(*
this);
604#endif
605
606 }
607
608
609 fitprob = (ncon - nunm > 0) ? gsl_cdf_chisq_Q(chinew, ncon - nunm) : 0.5;
610 chi2 = chinew;
611
612
613
614
615
616 gsl_matrix_set_zero(Vnew);
617 gsl_matrix_set_zero(Minv);
618 if (debug > 2) {
619 for (int i = 0; i < npar; ++i) {
620 for (int j = 0; j < npar; ++j) {
621 B2INFO("Minv[" << i << "," << j << "]=" << gsl_matrix_get(Minv, i, j));
622 }
623 }
624 }
625
626 B2DEBUG(10, "OPALFitterGSL: calcerr = " << calcerr);
627
628 if (calcerr) {
629
630
631
632
633
634
635
636 if (debug > 2) {
637 debug_print(&Vetaeta.matrix, "V");
638 debug_print(&Feta.matrix, "Feta");
639 }
640
641
642
643 gsl_blas_dsymm(CblasRight, CblasUpper, 1, &Vetaeta.matrix, &Feta.matrix, 0, FetaV);
644
645 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, FetaV, &Feta.matrix, 0, S);
646
647
648 if (nunm > 0) {
649
650
651
652
653
654 gsl_matrix_view Fxi = gsl_matrix_submatrix(Fetaxi, 0, nmea, ncon, nunm);
655
656 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, &Fxi.matrix, &Fxi.matrix, 1, S);
657 }
658
659 if (debug > 2) debug_print(S, "S");
660
661
662
663
664 int signum;
665 gsl_linalg_LU_decomp(S,
permS, &signum);
666 int inverr = gsl_linalg_LU_invert(S,
permS, Sinv);
667
668 if (inverr != 0) {
669 B2ERROR("S: gsl_linalg_LU_invert error " << inverr << " in error calculation");
670 ierr = -1;
671 return -1;
672 }
673
674
675
676
677
678
679
680 gsl_blas_dsymm(CblasLeft, CblasUpper, 1, Sinv, &Feta.matrix, 0, SinvFeta);
681
682 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1, &Feta.matrix, SinvFeta, 0, G);
683
684 if (debug > 2) debug_print(G, "G(1)");
685
686
687 if (nunm > 0) {
688
689
690
691
692 gsl_matrix_view Fxi = gsl_matrix_submatrix(Fetaxi, 0, nmea, ncon, nunm);
693
694
695 gsl_blas_dsymm(CblasLeft, CblasUpper, 1, Sinv, &Fxi.matrix, 0, SinvFxi);
696
697 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1, &Feta.matrix, SinvFxi, 0, H);
698
699 if (debug > 2) debug_print(H, "H");
700
701
702
703
704
705 gsl_matrix* Uinv = W1;
706 gsl_matrix_view U = gsl_matrix_submatrix(Minv, nmea, nmea, nunm, nunm);
707
708
709 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1, &Fxi.matrix, SinvFxi, 0, Uinv);
710
711 gsl_linalg_LU_decomp(Uinv,
permU, &signum);
712 inverr = gsl_linalg_LU_invert(Uinv,
permU, &U.matrix);
713
714 if (debug > 2) debug_print(&U.matrix, "U");
715 for (int i = 0; i < npar; ++i) {
716 for (int j = 0; j < npar; ++j) {
717 B2DEBUG(12, "after U Minv[" << i << "," << j << "]=" << gsl_matrix_get(Minv, i, j));
718 }
719 }
720
721 if (inverr != 0) {
722 B2ERROR("U: gsl_linalg_LU_invert error " << inverr << " in error calculation ");
723 ierr = -1;
724 return -1;
725 }
726
727
728
729
730 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, H, &U.matrix, 0, HU);
731
732 gsl_matrix_view Minvetaxi = gsl_matrix_submatrix(Minv, 0, nmea, nmea, nunm);
733 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, -1, &Vetaeta.matrix, HU, 0, &Minvetaxi.matrix);
734 for (int i = 0; i < npar; ++i) {
735 for (int j = 0; j < npar; ++j) {
736 B2DEBUG(12, "after etaxi Minv[" << i << "," << j << "]=" << gsl_matrix_get(Minv, i, j));
737 }
738 }
739
740
741
742
743 gsl_matrix_view Minvxieta = gsl_matrix_submatrix(Minv, nmea, 0, nunm, nmea);
744 gsl_matrix_transpose_memcpy(&Minvxieta.matrix, &Minvetaxi.matrix);
745 for (int i = 0; i < npar; ++i) {
746 for (int j = 0; j < npar; ++j) {
747 B2DEBUG(12, "after symmetric: Minv[" << i << "," << j << "]=" << gsl_matrix_get(Minv, i, j));
748 }
749 }
750
751
752
753 gsl_blas_dgemm(CblasNoTrans, CblasTrans, -1, HU, H, +1, G);
754
755 }
756
757
758
759 gsl_matrix_set_identity(IGV);
760
761 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, -1, G, &Vetaeta.matrix, 1, IGV);
762
763
764 gsl_matrix_view Minvetaeta = gsl_matrix_submatrix(Minv, 0, 0, nmea, nmea);
765
766
767 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &Vetaeta.matrix, IGV, 0, &Minvetaeta.matrix);
768
769 for (int i = 0; i < npar; ++i) {
770 for (int j = 0; j < npar; ++j) {
771 B2DEBUG(12, "complete Minv[" << i << "," << j << "]=" << gsl_matrix_get(Minv, i, j));
772 }
773 }
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791 gsl_matrix_view detadt = gsl_matrix_submatrix(dxdt, 0, 0, nmea, nmea);
792 B2DEBUG(13, "after detadt");
793
794 gsl_matrix_view Vdetadt = gsl_matrix_submatrix(Vdxdt, 0, 0, nmea, nmea);
795 B2DEBUG(13, "after Vdetadt");
796
797
798 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &Minvetaeta.matrix, Vinv, 0, &detadt.matrix);
799 if (debug > 2) debug_print(&detadt.matrix, "deta/dt");
800
801
802 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, &Vetaeta.matrix, &detadt.matrix, 0, &Vdetadt.matrix);
803 if (debug > 2) debug_print(&Vdetadt.matrix, "Vetata * deta/dt");
804
805 gsl_matrix_view Vnewetaeta = gsl_matrix_submatrix(Vnew, 0, 0, nmea, nmea);
806
807 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &detadt.matrix, &Vdetadt.matrix, 0, &Vnewetaeta.matrix);
808
809 if (debug > 2) debug_print(Vnew, "Vnew after part for measured parameters");
810
811
812 if (nunm > 0) {
813
814 gsl_matrix_view Minvxieta = gsl_matrix_submatrix(Minv, nmea, 0, nunm, nmea);
815 B2DEBUG(13, "after Minvxieta");
816 if (debug > 2) debug_print(&Minvxieta.matrix, "Minvxieta");
817
818 gsl_matrix_view dxidt = gsl_matrix_submatrix(dxdt, nmea, 0, nunm, nmea);
819 B2DEBUG(13, "after dxidt");
820
821 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &Minvxieta.matrix, Vinv, 0, &dxidt.matrix);
822 if (debug > 2) debug_print(&dxidt.matrix, "dxi/dt");
823
824
825 gsl_matrix_view Vdxidt = gsl_matrix_submatrix(Vdxdt, 0, nmea, nmea, nunm);
826 B2DEBUG(13, "after Vdxidt");
827
828 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1, &Vetaeta.matrix, &dxidt.matrix, 0, &Vdxidt.matrix);
829 if (debug > 2) debug_print(&Vdxidt.matrix, "Vetaeta * dxi/dt^T");
830
831 gsl_matrix_view Vnewetaxi = gsl_matrix_submatrix(Vnew, 0, nmea, nmea, nunm);
832 gsl_matrix_view Vnewxieta = gsl_matrix_submatrix(Vnew, nmea, 0, nunm, nmea);
833 gsl_matrix_view Vnewxixi = gsl_matrix_submatrix(Vnew, nmea, nmea, nunm, nunm);
834
835
836 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &dxidt.matrix, &Vdetadt.matrix, 0, &Vnewxieta.matrix);
837 if (debug > 2) debug_print(Vnew, "Vnew after xieta part");
838
839 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &detadt.matrix, &Vdxidt.matrix, 0, &Vnewetaxi.matrix);
840 if (debug > 2) debug_print(Vnew, "Vnew after etaxi part");
841
842 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1, &dxidt.matrix, &Vdxidt.matrix, 0, &Vnewxixi.matrix);
843 if (debug > 2) debug_print(Vnew, "Vnew after xixi part");
844 }
845
846
847
848
849 for (auto& fitobject : fitobjects) {
850 for (int ilocal = 0; ilocal < fitobject->getNPar(); ++ilocal) {
851 int iglobal = fitobject->getGlobalParNum(ilocal);
852 for (int jlocal = ilocal; jlocal < fitobject->getNPar(); ++jlocal) {
853 int jglobal = fitobject->getGlobalParNum(jlocal);
854 if (iglobal >= 0 && jglobal >= 0)
855 fitobject->setCov(ilocal, jlocal, gsl_matrix_get(Vnew, iglobal, jglobal));
856 }
857 }
858 }
859
860
864 }
867 for (
int i = 0; i <
covDim; ++i) {
868 for (
int j = 0; j <
covDim; ++j) {
869 cov[i *
covDim + j] = gsl_matrix_get(Vnew, i, j);
870 }
871 }
873
874 }
875
876#ifndef FIT_TRACEOFF
877 if (tracer) tracer->
finish(*
this);
878#endif
879
880 return fitprob;
881
882 }
double * cov
global covariance matrix of last fit problem
int covDim
dimension of global covariance matrix
virtual void step(BaseFitter &fitter)
Called at the end of each step.
virtual void initialize(BaseFitter &fitter)
Called at the start of a new fit (during initialization)
virtual void finish(BaseFitter &fitter)
Called at the end of a fit.