X-shooter Pipeline Reference Manual 3.8.15
xsh_utils_polynomial.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-Shooter Pipeline *
3 * Copyright (C) 2004,2005 European Southern Observatory *
4 * *
5 * This library is free software; you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation; either version 2 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2013-03-15 11:50:04 $
23 * $Revision: 1.8 $
24 * $Name: not supported by cvs2svn $
25 *
26 */
27
28#ifdef HAVE_CONFIG_H
29# include <config.h>
30#endif
31
32/*----------------------------------------------------------------------------*/
52/*----------------------------------------------------------------------------*/
53
54/*-----------------------------------------------------------------------------
55 Defines
56 -----------------------------------------------------------------------------*/
57
58/*
59 * When storing a 2d polynomial in a table,
60 * these column names are used
61 */
62#define COLUMN_ORDER1 "Order1"
63#define COLUMN_ORDER2 "Order2"
64#define COLUMN_COEFF "Coeff"
67/*-----------------------------------------------------------------------------
68 Includes
69 -----------------------------------------------------------------------------*/
71
72#include <xsh_utils.h>
73//#include <xsh_utils_wrappers.h>
74#include <xsh_dump.h>
75#include <xsh_msg.h>
76#include <xsh_error.h>
77
78#include <cpl.h>
79
80/*-----------------------------------------------------------------------------
81 Typedefs
82 -----------------------------------------------------------------------------*/
86{
88 cpl_polynomial *pol;
89
91 cpl_vector *vec;
92 double *vec_data;
93
94 int dimension; /* for efficiency */
95
97 double *shift;
98
100 double *scale;
101};
102
103
104
105/* COPY FROM CPL6.2 WITH sum changed from double to long double: DFS12436 */
106cpl_matrix * xsh_matrix_product_normal_create(const cpl_matrix * self)
107{
108
109 long double sum;
110 cpl_matrix * product;
111 const double * ai = cpl_matrix_get_data_const(self);
112 const double * aj;
113 double * bwrite;
114 const size_t m = cpl_matrix_get_nrow(self);
115 const size_t n = cpl_matrix_get_ncol(self);
116 size_t i, j, k;
117
118
119 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
120
121#if 0
122 /* Initialize all values to zero.
123 This is done to avoid access of uninitilized memory, in case
124 someone passes the matrix to for example cpl_matrix_dump(). */
125 product = cpl_matrix_new(m, m);
126 bwrite = cpl_matrix_get_data(product);
127#else
128 bwrite = (double *) cpl_malloc((size_t)m * (size_t)m * sizeof(double));
129 product = cpl_matrix_wrap(m, m, bwrite);
130#endif
131
132 /* The result at (i,j) is the dot-product of i'th and j'th row */
133 for (i = 0; i < m; i++, bwrite += m, ai += n) {
134 aj = ai; /* aj points to first entry in j'th row */
135 for (j = i; j < m; j++, aj += n) {
136 sum = 0.0;
137 for (k = 0; k < n; k++) {
138 sum += (long double)ai[k] * (long double)aj[k];
139 }
140 bwrite[j] = sum;
141 }
142 }
143
144 //cpl_tools_add_flops((cpl_flops)(n * m * (m + 1)));
145
146 return product;
147
148}
149
150/* COPY FROM CPL6.2 calling xsh_matrix_product_normal_create instead of cpl_matrix_product_normal_create
151 * and private symbols replaced with public ones */
152cpl_matrix *xsh_matrix_solve_normal(const cpl_matrix *coeff,
153 const cpl_matrix *rhs)
154{
155
156 cpl_matrix * solution;
157 cpl_matrix * At;
158 cpl_matrix * AtA;
159
160 cpl_ensure(coeff != NULL, CPL_ERROR_NULL_INPUT, NULL);
161 cpl_ensure(rhs != NULL, CPL_ERROR_NULL_INPUT, NULL);
162 //cpl_ensure(rhs->nr == coeff->nr, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
163
164 At = cpl_matrix_transpose_create(coeff);
165 solution = cpl_matrix_product_create(At, rhs);
166
168
169 cpl_matrix_delete(At);
170
171
172 /* not a public symbol, copy the two lines it does
173 if (cpl_matrix_solve_spd(AtA, solution)) {
174 */
175 if (cpl_matrix_decomp_chol(AtA) != CPL_ERROR_NONE ||
176 cpl_matrix_solve_chol(AtA, solution) != CPL_ERROR_NONE) {
177 cpl_matrix_delete(solution);
178 solution = NULL;
179 //(void)cpl_error_set_where_();
180 (void)cpl_error_set_where(cpl_func);
181 }
182
183 cpl_matrix_delete(AtA);
184
185 return solution;
186
187}
188
189
190
191/*-----------------------------------------------------------------------------
192 Implementation
193 -----------------------------------------------------------------------------*/
194/*----------------------------------------------------------------------------*/
205/*----------------------------------------------------------------------------*/
207xsh_polynomial_new(const cpl_polynomial *pol)
208{
209 polynomial *p = NULL;
210 int i;
211
212 /* Test input */
213 assure(pol != NULL, CPL_ERROR_ILLEGAL_INPUT, "Null polynomial");
214
215 /* Allocate and initialize struct */
216 p = cpl_calloc(1, sizeof(polynomial)) ;
217 assure_mem( p );
218
219 check_msg( p->dimension = cpl_polynomial_get_dimension(pol), "Error reading dimension");
220
221 /* Allocate vector */
222 p->vec = cpl_vector_new(p->dimension);
223 assure_mem( p->vec );
224 p->vec_data = cpl_vector_get_data(p->vec);
225
226 /* Shifts are initialized to zero, scales to 1 */
227 p->shift = cpl_calloc(p->dimension + 1, sizeof(double));
228 assure_mem( p->shift );
229
230 p->scale = cpl_malloc((p->dimension + 1) * sizeof(double));
231 assure_mem( p->scale );
232 for (i = 0; i <= p->dimension; i++)
233 p->scale[i] = 1.0;
234
235 check_msg( p->pol = cpl_polynomial_duplicate(pol), "Error copying polynomial");
236
237 cleanup:
238 if (cpl_error_get_code() != CPL_ERROR_NONE)
240
241 return p;
242}
243
244/*----------------------------------------------------------------------------*/
252/*----------------------------------------------------------------------------*/
255{
256 polynomial *result = NULL;
257 cpl_polynomial *p = NULL;
258
259 assure( dim >= 1, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dim);
260
261 p = cpl_polynomial_new(dim);
262 assure_mem( p );
263
264 result = xsh_polynomial_new(p);
265 assure_mem( result );
266
267 cleanup:
269
270 return result;
271}
272
273/*----------------------------------------------------------------------------*/
280/*----------------------------------------------------------------------------*/
281void
283{
285}
286
287/*----------------------------------------------------------------------------*/
294/*----------------------------------------------------------------------------*/
295void
297{
298 if (*p == NULL) return;
299 cpl_polynomial_delete((*p)->pol);
300 cpl_vector_delete((*p)->vec);
301 cpl_free((*p)->shift);
302 cpl_free((*p)->scale);
303 xsh_free(*p);
304 *p = NULL;
305 return;
306}
307/*----------------------------------------------------------------------------*/
313/*----------------------------------------------------------------------------*/
314int
316{
317 int result = -1;
318 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
319
320 result = cpl_polynomial_get_degree(p->pol);
321
322 cleanup:
323 return result;
324}
325
326/*----------------------------------------------------------------------------*/
332/*----------------------------------------------------------------------------*/
335{
336 polynomial *result = NULL;
337 int dimension;
338 int i;
339
340 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
342
343 check_msg( result = xsh_polynomial_new(p->pol),
344 "Error allocating polynomial");
345
346 for (i = 0; i <= dimension; i++)
347 {
348 result->shift[i] = p->shift[i];
349 result->scale[i] = p->scale[i];
350 }
351
352 cleanup:
353 if (cpl_error_get_code() != CPL_ERROR_NONE)
354 {
355 xsh_polynomial_delete(&result);
356 return NULL;
357 }
358
359 return result;
360}
361
362
363/*----------------------------------------------------------------------------*/
374/*----------------------------------------------------------------------------*/
375cpl_table *
377{
378 cpl_table *t = NULL; /* Result */
379 int degree;
380 int i, j, row;
381
382 /* Check input */
383 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
385 CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2D");
386
387 degree = cpl_polynomial_get_degree(p->pol);
388
389 /* Allocate space for 3 shifts, 3 scale factors and all
390 coefficients */
391 t = cpl_table_new(3 + 3 + (degree + 1)*(degree + 2)/2);
392 cpl_table_new_column(t, COLUMN_ORDER1, CPL_TYPE_INT);
393 cpl_table_new_column(t, COLUMN_ORDER2, CPL_TYPE_INT);
394 cpl_table_new_column(t, COLUMN_COEFF , CPL_TYPE_DOUBLE);
395
396 row = 0;
397
398 /* First write the shifts, write non-garbage to coeff columns (which are not used) */
399 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
400 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
401 cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[0]); row++;
402
403 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
404 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
405 cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[1]); row++;
406
407 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
408 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
409 cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[2]); row++;
410
411 /* Then the scale factors */
412 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
413 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
414 cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[0]); row++;
415
416 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
417 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
418 cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[1]); row++;
419
420 cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
421 cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
422 cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[2]); row++;
423
424 /* And then write the coefficients */
425 for (i = 0; i <= degree; i++){
426 for (j = 0; j+i <= degree; j++){
427 double coeff;
428 cpl_size power[2];
429 power[0] = i;
430 power[1] = j;
431
432 coeff = cpl_polynomial_get_coeff(p->pol, power);
433 cpl_table_set_int (t, COLUMN_ORDER1, row, power[0]);
434 cpl_table_set_int (t, COLUMN_ORDER2, row, power[1]);
435 cpl_table_set_double(t, COLUMN_COEFF , row, coeff);
436
437 row++;
438 }
439 }
440
441 cleanup:
442 return t;
443}
444
445/*----------------------------------------------------------------------------*/
454/*----------------------------------------------------------------------------*/
457{
458 polynomial *p = NULL; /* Result */
459 cpl_polynomial *pol = NULL;
460 cpl_type type;
461 int i;
462
463 /* Only 2d supported */
464 check_msg( pol = cpl_polynomial_new(2), "Error initializing polynomial");
465
466 /* Check table format */
467 assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
468 assure(cpl_table_has_column(t, COLUMN_ORDER1), CPL_ERROR_ILLEGAL_INPUT,
469 "No '%s' column found in table", COLUMN_ORDER1);
470 assure(cpl_table_has_column(t, COLUMN_ORDER2), CPL_ERROR_ILLEGAL_INPUT,
471 "No '%s' column found in table", COLUMN_ORDER2);
472 assure(cpl_table_has_column(t, COLUMN_COEFF ), CPL_ERROR_ILLEGAL_INPUT,
473 "No '%s' column found in table", COLUMN_COEFF );
474
475 type = cpl_table_get_column_type(t, COLUMN_ORDER1);
476 assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
477 "Column '%s' has type %s. Integer expected", COLUMN_ORDER1,
479
480 type = cpl_table_get_column_type(t, COLUMN_ORDER2);
481 assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
482 "Column '%s' has type %s. Integer expected", COLUMN_ORDER2,
484
485 type = cpl_table_get_column_type(t, COLUMN_COEFF);
486 assure(type == CPL_TYPE_DOUBLE, CPL_ERROR_INVALID_TYPE,
487 "Column '%s' has type %s. Double expected", COLUMN_COEFF ,
489
490 assure(cpl_table_get_nrow(t) > 1 + 2 + 1 + 2, CPL_ERROR_ILLEGAL_INPUT,
491 "Table must contain at least one coefficient");
492
493 /* Read the coefficients */
494 for(i = 3 + 3; i < cpl_table_get_nrow(t); i++) {
495 double coeff;
496 cpl_size power[2];
497
498 check_msg(( power[0] = cpl_table_get_int(t, COLUMN_ORDER1, i, NULL),
499 power[1] = cpl_table_get_int(t, COLUMN_ORDER2, i, NULL),
500 coeff = cpl_table_get_double(t, COLUMN_COEFF , i, NULL)),
501 "Error reading table row %d", i);
502
503 xsh_msg_debug("Pol.coeff.(%" CPL_SIZE_FORMAT ", %" CPL_SIZE_FORMAT ") = %e", power[0], power[1], coeff);
504
505 check_msg( cpl_polynomial_set_coeff(pol, power, coeff), "Error creating polynomial");
506 }
507 p = xsh_polynomial_new(pol);
508
509 /* Read shifts and rescaling */
510 xsh_polynomial_rescale(p, 0, cpl_table_get_double( t, COLUMN_COEFF, 3, NULL));
511 xsh_polynomial_rescale(p, 1, cpl_table_get_double( t, COLUMN_COEFF, 4, NULL));
512 xsh_polynomial_rescale(p, 2, cpl_table_get_double( t, COLUMN_COEFF, 5, NULL));
513 xsh_polynomial_shift (p, 0, cpl_table_get_double( t, COLUMN_COEFF, 0, NULL));
514 xsh_polynomial_shift (p, 1, cpl_table_get_double( t, COLUMN_COEFF, 1, NULL));
515 xsh_polynomial_shift (p, 2, cpl_table_get_double( t, COLUMN_COEFF, 2, NULL));
516
517 cleanup:
519 if (cpl_error_get_code() != CPL_ERROR_NONE)
521
522 return p;
523}
524
525
526/*----------------------------------------------------------------------------*/
532/*----------------------------------------------------------------------------*/
533int
535{
536 int dim = -1;
537 assure(p != NULL, CPL_ERROR_ILLEGAL_INPUT, "Null polynomial");
538
539/* slow check_msg( dim = cpl_polynomial_get_dimension(p->pol), "Error reading dimension"); */
540 dim = p->dimension;
541
542 cleanup:
543 return dim;
544}
545
546/*----------------------------------------------------------------------------*/
554/*----------------------------------------------------------------------------*/
555void xsh_polynomial_dump(const polynomial *p, FILE *stream)
556{
557 if (p == NULL)
558 fprintf(stream, "Null polynomial\n");
559 else {
560 int i;
561 cpl_polynomial_dump(p->pol, stream);
562 fprintf(stream, "shift_y \t= %f \tscale_y \t= %f\n", p->shift[0], p->scale[0]);
563 for (i = 1; i <= xsh_polynomial_get_dimension(p); i++)
564 {
565 fprintf(stream, "shift_x%d \t= %f \tscale_x%d \t= %f\n",
566 i, p->shift[i], i, p->scale[i]);
567 }
568 }
569 return;
570}
571
572/*----------------------------------------------------------------------------*/
586/*----------------------------------------------------------------------------*/
587cpl_error_code
588xsh_polynomial_rescale(polynomial *p, int varno, double scale)
589{
590 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
591 assure(0 <= varno && varno <= xsh_polynomial_get_dimension(p),
592 CPL_ERROR_ILLEGAL_INPUT, "Illegal variable number: %d", varno);
593
594 /* Rescaling an x variable by the factor S corresponds to:
595 * p'(x) := p(x/S) =
596 * cpl( (x/S - shiftx ) / scalex ) * scaley + shifty =
597 * cpl( (x - (S*shiftx)) / (S*scalex) ) * scaley + shifty */
598
599 /* Rescaling the y variable by the factor S corresponds to:
600 * p'(x) := S*p(x) =
601 * S * ( cpl((x - shiftx)/scalex) * scaley + shifty ) =
602 * cpl((x - shiftx)/scalex) * (S*scaley) + (S*shifty)
603 *
604 * therefore the implementation is the same in the two cases. */
605
606 p->shift[varno] *= scale;
607 p->scale[varno] *= scale;
608
609 cleanup:
610 return cpl_error_get_code();
611}
612
613/*----------------------------------------------------------------------------*/
627/*----------------------------------------------------------------------------*/
628cpl_error_code
629xsh_polynomial_shift(polynomial *p, int varno, double shift)
630{
631 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
632 assure(0 <= varno && varno <= xsh_polynomial_get_dimension(p),
633 CPL_ERROR_ILLEGAL_INPUT, "Illegal variable number: %d", varno);
634
635 /* The implementation is similar for x and y variables because
636 * p(x-S) =
637 * cpl( (x-S - shiftx) / scalex ) * scaley + shifty =
638 * cpl( (x - (shiftx+S)) / scalex ) * scaley + shifty
639 * and
640 * p(x) + S =
641 * cpl( (x - shiftx)/scalex ) * scaley + shifty + S =
642 * cpl( (x - shiftx)/scalex ) * scaley + (shifty+S) */
643
644 p->shift[varno] += shift;
645
646 cleanup:
647 return cpl_error_get_code();
648}
649
650/*----------------------------------------------------------------------------*/
659/*----------------------------------------------------------------------------*/
660double
662{
663 double result = 0;
664
665 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
667 CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 1d");
668
669 check_msg( result =
670 cpl_polynomial_eval_1d(p->pol, (x - p->shift[1])/p->scale[1], NULL)
671 * p->scale[0] + p->shift[0],
672 "Could not evaluate polynomial");
673
674 cleanup:
675 return result;
676}
677
678
679/*----------------------------------------------------------------------------*/
689/*----------------------------------------------------------------------------*/
690
691double
692xsh_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
693{
694 double result = 0;
695
696 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
697 assure(p->dimension == 2, CPL_ERROR_ILLEGAL_INPUT,
698 "Polynomial must be 2d. It's %dd", p->dimension);
699 {
700 double scale = p->scale[0];
701 double shift = p->shift[0];
702
703 // cpl_vector_set(p->vec, 0, (x1 - p->shift[1]) / p->scale[1]);
704 // cpl_vector_set(p->vec, 1, (x2 - p->shift[2]) / p->scale[2]);
705 p->vec_data[0] = (x1 - p->shift[1]) / p->scale[1];
706 p->vec_data[1] = (x2 - p->shift[2]) / p->scale[2];
707
708 result = cpl_polynomial_eval(p->pol, p->vec) * scale + shift;
709 }
710
711 cleanup:
712 return result;
713}
714
715/*----------------------------------------------------------------------------*/
728/*----------------------------------------------------------------------------*/
729double
730xsh_polynomial_solve_1d(const polynomial *p, double value, double guess, int multiplicity)
731{
732 double result = 0;
733 cpl_size power[1];
734 double coeff0;
735
736 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
737 assure(xsh_polynomial_get_dimension(p) == 1, CPL_ERROR_ILLEGAL_INPUT,
738 "Polynomial must be 1d");
739
740 /* Solving p(x) = value corresponds to solving
741 <=> cpl_p( (x-xshift)/xscale )*yscale + yshift = value
742 <=> cpl_p( (x-xshift)/xscale ) + (yshift - value)/yscale = 0
743
744 So 1) find zero point for the polynomial cpl() + (yshift-value)/yscale
745 Then 2) shift and rescale the result
746 */
747
748 power[0] = 0;
749 check_msg(( coeff0 = cpl_polynomial_get_coeff(p->pol, power),
750 cpl_polynomial_set_coeff(p->pol, power, coeff0 + (p->shift[0] - value)/p->scale[0])),
751 "Error setting coefficient");
752
753 check_msg( cpl_polynomial_solve_1d(p->pol, (guess - p->shift[1]) / p->scale[1],
754 &result, multiplicity), "Could not find root");
755 /* Restore polynomial */
756 cpl_polynomial_set_coeff(p->pol, power, coeff0);
757
758 /* Shift solution */
759 result = result * p->scale[1] + p->shift[1];
760
761 cleanup:
762 return result;
763}
764
765/*----------------------------------------------------------------------------*/
782/*----------------------------------------------------------------------------*/
783double
784xsh_polynomial_solve_2d(const polynomial *p, double value, double guess,
785 int multiplicity, int varno, double x_value)
786{
787 double result = 0;
788 polynomial *pol_1d = NULL;
789
790 assure( 1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
791 "Illegal variable number: %d", varno);
792
793 check_msg( pol_1d = xsh_polynomial_collapse(p, varno, x_value),
794 "Could not collapse polynomial");
795
796 check_msg( result = xsh_polynomial_solve_1d(pol_1d, value, guess, multiplicity),
797 "Could not find root");
798
799 cleanup:
800 xsh_polynomial_delete(&pol_1d);
801 return result;
802}
803
804/*----------------------------------------------------------------------------*/
813/*----------------------------------------------------------------------------*/
814double
815xsh_polynomial_derivative_2d(const polynomial *p, double x1, double x2, int varno)
816{
817 double result = 0;
818 cpl_size power[2];
819
820 assure (1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
821 "Illegal variable number (%d)", varno);
822
823 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
824 assure(xsh_polynomial_get_dimension(p) == 2, CPL_ERROR_ILLEGAL_INPUT,
825 "Polynomial must be 2d. It's %dd", xsh_polynomial_get_dimension(p));
826
827 /* d/dx_i [ p(x) ] =
828 * d/dx_i [ cpl( (x - shiftx) / scalex ) * scaley + shifty ] =
829 * [ d(cpl)/dx_i ( (x - shiftx) / scalex ) * scaley ]
830 */
831
832 /* Shift, scale (x1, x2) */
833 x1 = (x1 - p->shift[1])/p->scale[1];
834 x2 = (x2 - p->shift[2])/p->scale[2];
835
836 /* Get derivative of cpl polynomial.
837 *
838 */
839 {
840 int degree = cpl_polynomial_get_degree(p->pol);
841 double yj = 1; /* y^j */
842 int i, j;
843
844 result = 0;
845 for (j = 0, yj = 1;
846 j <= degree; j++,
847 yj *= (varno == 1) ? x2 : x1)
848 {
849 /* Proof by example (degree = 3): For each j account for these terms
850 * using Horner's rule:
851 *
852 * d/dx y^j * [ c_3j x^3 + c_2j x^2 + c_1j x^1 + c_0j ] =
853 *
854 * y^j * [ 3c_3j x^2 + 2c_2j x^1 + 1c_1j ] =
855 *
856 * y^j * [ ((3c_3j) x + 2c_2j) x + 1c_1j ]
857 */
858
859 double sum = 0;
860 for (i = degree; i >= 1; i--)
861 {
862 double c_ij;
863
864 power[0] = (varno == 1) ? i : j;
865 power[1] = (varno == 1) ? j : i;
866
867 c_ij = cpl_polynomial_get_coeff(p->pol, power);
868
869 sum += (i * c_ij);
870 if (i >= 2) sum *= (varno == 1) ? x1 : x2;
871 }
872
873 /* Collect terms */
874 result += yj * sum;
875 }
876 }
877
878 result *= p->scale[0];
879
880
881/* Old code: This method (valid for varno = 2)
882 of getting the derivative of
883 the CPL polynomial is slow because of the call to
884 xsh_polynomial_collapse()
885
886 check_msg( pol_1d = xsh_polynomial_collapse(p, 1, x1);
887 dummy = cpl_polynomial_eval_1d(pol_1d->pol, (x2 - pol_1d->shift[1])/pol_1d->scale[1], &result),
888 "Error evaluating derivative");
889*/
890
891 cleanup:
892 return result;
893}
894
895/*----------------------------------------------------------------------------*/
902/*----------------------------------------------------------------------------*/
903double
905{
906 double result = 0;
907
908 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
910 CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 1d");
911
912 check_msg( cpl_polynomial_eval_1d(p->pol, (x - p->shift[1])/p->scale[1], &result),
913 "Error evaluating derivative");
914
915 cleanup:
916 return result;
917}
918
919/*----------------------------------------------------------------------------*/
926/*----------------------------------------------------------------------------*/
929{
930 polynomial *result = NULL;
931 cpl_polynomial *pol = NULL;
932
933 assure(p1 != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
934 assure(p2 != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
936 CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
938 CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
939
940 /* cpl_polynomial1((x - shift_x1)/scale_x1) * scale_y1 + shift_y1
941 +
942 cpl_polynomial2((x - shift_x2)/scale_x2) * scale_y2 + shift_y2
943 = ???
944 Not easy.
945
946 Use brute force:
947 */
948
949 {
950 int degree, i, j;
951
954
955 pol = cpl_polynomial_new(2);
956 for (i = 0; i <= degree; i++)
957 for (j = 0; j <= degree; j++) {
958 double coeff1, coeff2;
959 cpl_size power[2];
960
961 /* Simple: add coefficients of the same power */
962 coeff1 = xsh_polynomial_get_coeff_2d(p1, i, j);
963 coeff2 = xsh_polynomial_get_coeff_2d(p2, i, j);
964
965 power[0] = i;
966 power[1] = j;
967 cpl_polynomial_set_coeff(pol, power, coeff1 + coeff2);
968 }
969 }
970
971 result = xsh_polynomial_new(pol);
972
973 cleanup:
975 return result;
976}
977
978/*----------------------------------------------------------------------------*/
991/*----------------------------------------------------------------------------*/
992static cpl_error_code
993derivative_cpl_polynomial(cpl_polynomial *p, int varno)
994{
995 int dimension, degree;
996 int i, j;
997 cpl_size power[2];
998
999 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1000 dimension = cpl_polynomial_get_dimension(p);
1001 degree = cpl_polynomial_get_degree(p);
1002 assure( 1 <= dimension && dimension <= 2, CPL_ERROR_ILLEGAL_INPUT,
1003 "Illegal dimension: %d", dimension);
1004 assure( 1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1005 "Illegal variable number: %d", varno);
1006
1007 if (dimension == 1)
1008 {
1009 /* a_i := (i+1) * a_(i+1) */
1010 for(i = 0; i <= degree; i++)
1011 {
1012 double coeff;
1013 power[0] = i+1;
1014 /* power[1] is ignored */
1015
1016 coeff = cpl_polynomial_get_coeff(p, power);
1017
1018 power[0] = i;
1019 cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1020 }
1021 }
1022
1023 if (dimension == 2)
1024 {
1025 /* a_ij := (i+1) * a_{(i+1),j} */
1026 for(i = 0; i <= degree; i++)
1027 {
1028 for(j = 0; i + j <= degree; j++)
1029 {
1030 double coeff;
1031 power[varno - 1] = i+1; /* varno == 1: 0,1 */
1032 power[2 - varno] = j; /* varno == 2: 1,0 */
1033
1034 coeff = cpl_polynomial_get_coeff(p, power);
1035
1036 power[varno - 1] = i;
1037
1038 cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1039 }
1040 }
1041 }
1042
1043 cleanup:
1044 return cpl_error_get_code();
1045}
1046
1047/*----------------------------------------------------------------------------*/
1057/*----------------------------------------------------------------------------*/
1058cpl_error_code
1060{
1061 int dimension;
1062
1063 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1064 check_msg ( dimension = xsh_polynomial_get_dimension(p), "Error reading dimension");
1065 assure( 1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1066 "Illegal variable number: %d", varno);
1067
1068
1069 /* d/dx_i [ cpl( (x - shiftx) / scalex ) * scaley + shifty ] =
1070 * sum_j d(cpl)/dx_j ( (x - shiftx) / scalex ) * scaley * dx_j/dx_i / scalex_j =
1071 * d(cpl)/dx_i ( (x - shiftx) / scalex ) * scaley/scalex_i,
1072 *
1073 * so transform : shifty -> 0
1074 * shiftx -> shiftx
1075 * scaley -> scaley/scalex_i
1076 * scalex -> scalex
1077 * cpl -> d(cpl)/dx_i
1078 */
1079
1080 p->shift[0] = 0;
1081 p->scale[0] = p->scale[0] / p->scale[varno];
1082
1084 "Error calculating derivative of CPL-polynomial");
1085
1086 cleanup:
1087 return cpl_error_get_code();
1088}
1089
1090
1091/*----------------------------------------------------------------------------*/
1100/*----------------------------------------------------------------------------*/
1101double
1102xsh_polynomial_get_coeff_2d(const polynomial *p, int degree1, int degree2)
1103{
1104 polynomial *pp = NULL;
1105 int dimension;
1106 double result = 0;
1107 double factorial;
1108
1109 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1110 check_msg ( dimension = xsh_polynomial_get_dimension(p), "Error reading dimension");
1111 assure(dimension == 2, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dimension);
1112 assure( 0 <= degree1, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree1);
1113 assure( 0 <= degree2, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree2);
1114
1115 /* Calculate the coefficient as
1116 * d^N p / (dx1^degree1 dx2^degree2) / (degree1! * degree2!)
1117 * evaluated in (0,0)
1118 */
1119
1121
1122 factorial = 1;
1123 while(degree1 > 0)
1124 {
1125 check_msg( xsh_polynomial_derivative(pp, 1), "Error calculating derivative");
1126
1127 factorial *= degree1;
1128 degree1 -= 1;
1129 }
1130
1131 while(degree2 > 0)
1132 {
1133 check_msg( xsh_polynomial_derivative(pp, 2), "Error calculating derivative");
1134
1135 factorial *= degree2;
1136 degree2 -= 1;
1137 }
1138
1139 check_msg( result = xsh_polynomial_evaluate_2d(pp, 0, 0) / factorial,
1140 "Error evaluating polynomial");
1141
1142 cleanup:
1144 return result;
1145}
1146/*----------------------------------------------------------------------------*/
1156/*----------------------------------------------------------------------------*/
1157double
1159{
1160 polynomial *pp = NULL;
1161 int dimension;
1162 double result = 0;
1163 double factorial;
1164
1165 assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1166 check_msg ( dimension = xsh_polynomial_get_dimension(p), "Error reading dimension");
1167 assure(dimension == 1, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dimension);
1168 assure( 0 <= degree, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree);
1169
1170 /* Calculate the coefficient as
1171 * d^degree p/dx^degree / (degree1! * degree2!)
1172 * evaluated in 0.
1173 */
1174
1176
1177 factorial = 1;
1178 while(degree > 0)
1179 {
1180 check_msg( xsh_polynomial_derivative(pp, 1), "Error calculating derivative");
1181
1182 factorial *= degree;
1183 degree -= 1;
1184 }
1185
1186 check_msg( result = xsh_polynomial_evaluate_1d(pp, 0) / factorial,
1187 "Error evaluating polynomial");
1188
1189 cleanup:
1191 return result;
1192}
1193
1194
1195/*----------------------------------------------------------------------------*/
1211/*----------------------------------------------------------------------------*/
1212polynomial *
1213xsh_polynomial_collapse(const polynomial *p, int varno, double value)
1214{
1215 polynomial *result = NULL;
1216 cpl_polynomial *pol = NULL;
1217 cpl_size *power = NULL;
1218
1219 int i, j;
1220 int degree, dimension;
1221
1222 assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1224 assure(dimension > 0, CPL_ERROR_ILLEGAL_INPUT,
1225 "Polynomial has non-positive dimension: %d", dimension);
1226 assure(dimension != 1, CPL_ERROR_ILLEGAL_OUTPUT,
1227 "Don't collapse a 1d polynomial. Evaluate it!");
1228
1229 /* To gecpl_image_get_maxpos_windowneralize this function to work with dimensions higher than 2,
1230 also changes needs to be made below (use varno properly). For now,
1231 support only 2d. */
1232 assure(dimension == 2, CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
1233
1234 assure(1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1235 "Wrong variable number");
1236 value = (value - p->shift[varno]) / p->scale[varno];
1237
1238 /* Compute new coefficients */
1239 degree = cpl_polynomial_get_degree(p->pol);
1240 pol = cpl_polynomial_new(dimension - 1);
1241 power = cpl_malloc(sizeof(cpl_size) * dimension);
1242 assure_mem( power );
1243 for (i = 0; i <= degree; i++)
1244 {
1245 double coeff;
1246
1247 power[2-varno] = i; /* map 2->0 and 1->1 */
1248
1249 /* Collect all terms with x^i (using Horner's rule) */
1250 coeff = 0;
1251 for (j = degree - i; j >= 0; j--)
1252 {
1253 power[varno-1] = j; /* map 2->1 and 1->0 */
1254 coeff += cpl_polynomial_get_coeff(p->pol, power);
1255 if (j > 0) coeff *= value;
1256 }
1257 /* Write coefficient in 1d polynomial */
1258 power[0] = i;
1259 cpl_polynomial_set_coeff(pol, power, coeff);
1260 }
1261
1262 /* Wrap the polynomial */
1263 result = xsh_polynomial_new(pol);
1264
1265 /* Copy the shifts and scales, skip variable number varno */
1266 j = 0;
1267 for(i = 0; i <= dimension - 1; i++)
1268 {
1269 if (i == varno)
1270 {
1271 /* Don't copy */
1272 j += 2;
1273 /* For the remainder of this for loop, j = i+1 */
1274 }
1275 else
1276 {
1277 result->shift[i] = p->shift[j];
1278 result->scale[i] = p->scale[j];
1279 j += 1;
1280 }
1281 }
1282
1283 assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
1284 "Error collapsing polynomial");
1285
1286 cleanup:
1287 cpl_free(power); power = NULL;
1288 xsh_free_polynomial(&pol);
1289 if (cpl_error_get_code() != CPL_ERROR_NONE)
1290 {
1291 xsh_polynomial_delete(&result);
1292 }
1293 return result;
1294}
1295
1296
1297
1298/*----------------------------------------------------------------------------*/
1318/*----------------------------------------------------------------------------*/
1320 const cpl_vector * x_pos,
1321 const cpl_vector * values,
1322 const cpl_vector * sigmas,
1323 int poly_deg,
1324 double * mse)
1325{
1326 int nc ;
1327 int np ;
1328 cpl_matrix * ma = NULL;
1329 cpl_matrix * mb = NULL;
1330 cpl_matrix * mx = NULL;
1331 const double * x_pos_data ;
1332 const double * values_data ;
1333 const double * sigmas_data = NULL;
1334 double mean_x, mean_z;
1335 polynomial * result = NULL;
1336 cpl_polynomial * out ;
1337 cpl_vector * x_val = NULL;
1338 int i, j ;
1339
1340 /* Check entries */
1341 assure_nomsg( x_pos != NULL && values != NULL, CPL_ERROR_NULL_INPUT);
1342 assure( poly_deg >= 0, CPL_ERROR_ILLEGAL_INPUT,
1343 "Polynomial degree is %d. Must be non-negative", poly_deg);
1344 np = cpl_vector_get_size(x_pos) ;
1345
1346 nc = 1 + poly_deg ;
1347 assure( np >= nc, CPL_ERROR_ILLEGAL_INPUT,
1348 "Not enough points (%d) to fit %d-order polynomial. %d point(s) needed",
1349 np, poly_deg, nc);
1350
1351 /* Fill up look-up table for coefficients to compute */
1352 /* Initialize matrices */
1353 /* ma contains the polynomial terms for each input point. */
1354 /* mb contains the values */
1355 ma = cpl_matrix_new(np, nc) ;
1356 mb = cpl_matrix_new(np, 1) ;
1357
1358 /* Get mean values */
1359 mean_x = cpl_vector_get_mean(x_pos);
1360 mean_z = cpl_vector_get_mean(values);
1361
1362 /* Fill up matrices, shift */
1363 x_pos_data = cpl_vector_get_data_const(x_pos) ;
1364 values_data = cpl_vector_get_data_const(values) ;
1365 if (sigmas != NULL)
1366 {
1367 sigmas_data = cpl_vector_get_data_const(sigmas) ;
1368 }
1369
1370 if (sigmas != NULL)
1371 {
1372 for (i=0 ; i<np ; i++)
1373 {
1374 /* Catch division by zero */
1375 if (sigmas_data[i] == 0)
1376 {
1377 xsh_free_matrix(&ma) ;
1378 xsh_free_matrix(&mb) ;
1379 assure(false, CPL_ERROR_DIVISION_BY_ZERO,
1380 "Sigmas must be non-zero");
1381 }
1382 for (j=0 ; j<nc ; j++)
1383 {
1384 cpl_matrix_set(ma, i, j,
1385 xsh_pow_int(x_pos_data[i] - mean_x, j) /
1386 sigmas_data[i]) ;
1387 }
1388 /* mb contains surface values (z-axis) */
1389 cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / sigmas_data[i]);
1390 }
1391 }
1392 else /* Use sigma = 1 */
1393 {
1394 for (i=0 ; i<np ; i++)
1395 {
1396 for (j=0 ; j<nc ; j++)
1397 {
1398 cpl_matrix_set(ma, i, j,
1399 xsh_pow_int(x_pos_data[i] - mean_x, j) / 1);
1400 }
1401 /* mb contains surface values (z-values) */
1402 cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / 1) ;
1403 }
1404 }
1405
1406 /* Solve XA=B by a least-square solution (aka pseudo-inverse). */
1407 check_msg( mx = xsh_matrix_solve_normal(ma, mb),
1408 "Could not invert matrix");
1409 xsh_free_matrix(&ma);
1410 xsh_free_matrix(&mb);
1411
1412 /* Store coefficients for output */
1413 out = cpl_polynomial_new(1) ;
1414 cpl_size deg=0;
1415 for (deg=0 ; deg<nc ; deg++) {
1416 cpl_polynomial_set_coeff(out, &deg, cpl_matrix_get(mx, deg, 0)) ;
1417 }
1418 xsh_free_matrix(&mx);
1419
1420 /* If requested, compute mean squared error */
1421 if (mse != NULL) {
1422 *mse = 0.00 ;
1423 x_val = cpl_vector_new(1) ;
1424 for (i=0 ; i<np ; i++)
1425 {
1426 double residual;
1427 cpl_vector_set(x_val, 0, x_pos_data[i] - mean_x) ;
1428 /* Subtract from the true value, square, accumulate */
1429 residual = (values_data[i] - mean_z) - cpl_polynomial_eval(out, x_val);
1430 *mse += residual*residual;
1431 }
1432 xsh_free_vector(&x_val) ;
1433 /* Average the error term */
1434 *mse /= (double)np ;
1435 }
1436
1437 /* Create and shift result */
1438 result = xsh_polynomial_new(out);
1439 xsh_free_polynomial(&out);
1440
1441 xsh_polynomial_shift(result, 0, mean_z);
1442 xsh_polynomial_shift(result, 1, mean_x);
1443
1444 cleanup:
1445 xsh_free_vector(&x_val);
1446 xsh_free_matrix(&ma);
1447 xsh_free_matrix(&mb);
1448 xsh_free_matrix(&mx);
1449 return result;
1450}
1451
1452
1453/*----------------------------------------------------------------------------*/
1497/*----------------------------------------------------------------------------*/
1498polynomial *
1499xsh_polynomial_fit_2d(const cpl_bivector * xy_pos, const cpl_vector * values,
1500 const cpl_vector * sigmas, int poly_deg1, int poly_deg2, double * mse,
1501 double * red_chisq, polynomial ** variance) {
1502 int nc;
1503 int degx, degy;
1504 int * degx_tab;
1505 int * degy_tab;
1506 int np;
1507 cpl_matrix * ma;
1508 cpl_matrix * mb;
1509 cpl_matrix * mx;
1510 cpl_matrix * mat;
1511 cpl_matrix * mat_ma;
1512 cpl_matrix * cov = NULL;
1513 const double * xy_pos_data_x;
1514 const double * xy_pos_data_y;
1515 const double * values_data;
1516 const double * sigmas_data = NULL;
1517 const cpl_vector* xy_pos_x;
1518 const cpl_vector* xy_pos_y;
1519 double mean_x, mean_y, mean_z;
1520 cpl_polynomial * out;
1521 cpl_polynomial * variance_cpl;
1522 polynomial * result = NULL;
1523 cpl_size * powers;
1524 int i_nc=0;
1525
1526 /* Check entries */
1527 assure(xy_pos && values, CPL_ERROR_NULL_INPUT, "Null input");
1528 assure(poly_deg1 >= 0, CPL_ERROR_ILLEGAL_INPUT,
1529 "Polynomial degree1 is %d", poly_deg1);
1530 assure(poly_deg2 >= 0, CPL_ERROR_ILLEGAL_INPUT,
1531 "Polynomial degree2 is %d", poly_deg2);
1532 np = cpl_bivector_get_size(xy_pos);
1533
1534 /* Can't calculate variance and chi_sq without sigmas */
1535 assure( (variance == NULL && red_chisq == NULL) || sigmas != NULL,
1536 CPL_ERROR_ILLEGAL_INPUT,
1537 "Cannot calculate variance or chi_sq without knowing");
1538
1539 /* Fill up look-up table for coefficients to compute */
1540 nc = (1 + poly_deg1) * (1 + poly_deg2); /* rectangular matrix */
1541
1542 assure(np >= nc, CPL_ERROR_SINGULAR_MATRIX,
1543 "%d coefficients. Only %d points", nc, np);
1544 /* The error code here is set to SINGULAR_MATRIX, in order to allow the caller
1545 to detect when too many coefficients are fitted to too few points */
1546
1547 /* Need an extra point to calculate reduced chi^2 */
1548 assure(red_chisq == NULL || np > nc, CPL_ERROR_ILLEGAL_INPUT,
1549 "%d coefficients. %d points. Cannot calculate chi square", nc, np);
1550
1551 degx_tab = cpl_malloc(nc * sizeof(int));
1552 assure_mem( degx_tab);
1553
1554 degy_tab = cpl_malloc(nc * sizeof(int));
1555 if (degy_tab == NULL) {
1556 cpl_free(degx_tab);
1557 assure_mem( false);
1558 }
1559
1560 {
1561 int i = 0;
1562 for (degy = 0; degy <= poly_deg2; degy++) { /* rectangular matrix */
1563 for (degx = 0; degx <= poly_deg1; degx++) {
1564 degx_tab[i] = degx;
1565 degy_tab[i] = degy;
1566 i++;
1567 }
1568 }
1569 }
1570
1571 /* Initialize matrices */
1572 /* ma contains the polynomial terms in the order described */
1573 /* above in each column, for each input point. */
1574 /* mb contains the values */
1575 ma = cpl_matrix_new(np, nc);
1576 mb = cpl_matrix_new(np, 1);
1577
1578 /* Get the mean of each variable */
1579 xy_pos_x = cpl_bivector_get_x_const(xy_pos);
1580 xy_pos_y = cpl_bivector_get_y_const(xy_pos);
1581
1582 mean_x = cpl_vector_get_mean(xy_pos_x);
1583 mean_y = cpl_vector_get_mean(xy_pos_y);
1584 mean_z = cpl_vector_get_mean(values);
1585
1586 /* Fill up matrices. At the same time shift the data
1587 so that it is centered around zero */
1588 xy_pos_data_x = cpl_vector_get_data_const(xy_pos_x);
1589 xy_pos_data_y = cpl_vector_get_data_const(xy_pos_y);
1590 values_data = cpl_vector_get_data_const(values);
1591 if (sigmas != NULL) {
1592 sigmas_data = cpl_vector_get_data_const(sigmas);
1593 }
1594
1595 if (sigmas != NULL) {
1596 int i;
1597 for (i = 0; i < np; i++) {
1598 double *ma_data = cpl_matrix_get_data(ma);
1599 double *mb_data = cpl_matrix_get_data(mb);
1600
1601 int j = 0;
1602 double valy = 1;
1603
1604 /* Catch division by zero */
1605 if (sigmas_data[i] == 0) {
1606 xsh_free_matrix(&ma);
1607 xsh_free_matrix(&mb);
1608 cpl_free(degx_tab);
1609 cpl_free(degy_tab);
1610 assure(false, CPL_ERROR_DIVISION_BY_ZERO,
1611 "Sigmas must be non-zero. sigma[%d] is %f", i, sigmas_data[i]);
1612 }
1613 i_nc=i*nc;
1614 for (degy = 0; degy <= poly_deg2; degy++) {
1615 double valx = 1;
1616 for (degx = 0; degx <= poly_deg1; degx++) {
1617 ma_data[j + i_nc] = valx * valy / sigmas_data[i];
1618 valx *= (xy_pos_data_x[i] - mean_x);
1619 j++;
1620 }
1621 valy *= (xy_pos_data_y[i] - mean_y);
1622 }
1623
1624 /* mb contains surface values (z-axis) */
1625
1626 mb_data[i] = (values_data[i] - mean_z) / sigmas_data[i];
1627 }
1628 } else /* Use sigma = 1 */
1629 {
1630 int i;
1631 for (i = 0; i < np; i++) {
1632 double *ma_data = cpl_matrix_get_data(ma);
1633 double *mb_data = cpl_matrix_get_data(mb);
1634
1635 double valy = 1;
1636 i_nc=i*nc;
1637 int j = 0;
1638 for (degy = 0; degy <= poly_deg2; degy++) {
1639 double valx = 1;
1640 for (degx = 0; degx <= poly_deg1; degx++) {
1641 ma_data[j + i_nc] = valx * valy / 1;
1642 valx *= (xy_pos_data_x[i] - mean_x);
1643 j++;
1644 }
1645 valy *= (xy_pos_data_y[i] - mean_y);
1646 }
1647
1648 /* mb contains surface values (z-axis) */
1649// cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / 1) ;
1650 mb_data[i] = values_data[i] - mean_z;
1651 }
1652 }
1653
1654 /* If variance polynomial is requested,
1655 compute covariance matrix = (A^T * A)^-1 */
1656 if (variance != NULL) {
1657 mat = cpl_matrix_transpose_create(ma);
1658 if (mat != NULL) {
1659 mat_ma = cpl_matrix_product_create(mat, ma);
1660 if (mat_ma != NULL) {
1661 cov = cpl_matrix_invert_create(mat_ma);
1662 /* Here, one might do a (paranoia) check that the covariance
1663 matrix is symmetrical and has positive eigenvalues (so that
1664 the returned variance polynomial is guaranteed to be positive) */
1665
1666 variance_cpl = cpl_polynomial_new(2);
1667 }
1668 }
1670 xsh_free_matrix(&mat_ma);
1671 }
1672
1673 /* Solve XA=B by a least-square solution (aka pseudo-inverse). */
1674 mx = xsh_matrix_solve_normal(ma, mb);
1675
1676 xsh_free_matrix(&ma);
1677 xsh_free_matrix(&mb);
1678 if (mx == NULL) {
1679 cpl_free(degx_tab);
1680 cpl_free(degy_tab);
1681 xsh_free_matrix(&cov);
1682 assure(false, CPL_ERROR_ILLEGAL_OUTPUT, "Matrix inversion failed");
1683 }
1684
1685 /* Store coefficients for output */
1686 out = cpl_polynomial_new(2);
1687 powers = cpl_malloc(2 * sizeof(cpl_size));
1688 if (powers == NULL) {
1689 cpl_free(degx_tab);
1690 cpl_free(degy_tab);
1691 xsh_free_matrix(&mx);
1692 xsh_free_matrix(&cov);
1693 xsh_free_polynomial(&out);
1694 assure_mem( false);
1695 }
1696
1697 {
1698 int i;
1699 for (i = 0; i < nc; i++) {
1700 powers[0] = degx_tab[i];
1701 powers[1] = degy_tab[i];
1702 cpl_polynomial_set_coeff(out, powers, cpl_matrix_get(mx, i, 0));
1703
1704 /* Create variance polynomial (if requested) */
1705 if (variance != NULL && /* Requested? */
1706 cov != NULL && variance_cpl != NULL /* covariance computation succeeded? */
1707 ) {
1708 int j;
1709 for (j = 0; j < nc; j++) {
1710 double coeff;
1711 /* Add cov_ij to the proper coeff:
1712 cov_ij * dp/d(ai) * dp/d(aj) =
1713 cov_ij * (x^degx[i] * y^degy[i]) * (x^degx[i] * y^degy[i]) =
1714 cov_ij * x^(degx[i]+degx[j]) * y^(degy[i] + degy[j]),
1715
1716 i.e. add cov_ij to coeff (degx[i]+degx[j], degy[i]+degy[j]) */
1717 powers[0] = degx_tab[i] + degx_tab[j];
1718 powers[1] = degy_tab[i] + degy_tab[j];
1719
1720 coeff = cpl_polynomial_get_coeff(variance_cpl, powers);
1721 cpl_polynomial_set_coeff(variance_cpl, powers,
1722 coeff + cpl_matrix_get(cov, i, j));
1723 }
1724 }
1725 }
1726 }
1727
1728 cpl_free(powers);
1729 cpl_free(degx_tab);
1730 cpl_free(degy_tab);
1731 xsh_free_matrix(&cov);
1732 xsh_free_matrix(&mx);
1733
1734 /* Create and shift result */
1735 result = xsh_polynomial_new(out);
1736 xsh_free_polynomial(&out);
1737 xsh_polynomial_shift(result, 0, mean_z);
1738 xsh_polynomial_shift(result, 1, mean_x);
1739 xsh_polynomial_shift(result, 2, mean_y);
1740
1741 /* Wrap up variance polynomial */
1742 if (variance != NULL) {
1743 *variance = xsh_polynomial_new(variance_cpl);
1744 xsh_free_polynomial(&variance_cpl);
1745 /* The variance of the fit does not change
1746 when a constant is added to the a_00
1747 coefficient of the polynomial, so don't:
1748 xsh_polynomial_shift(*variance, 0, mean_z); */
1749 xsh_polynomial_shift(*variance, 1, mean_x);
1750 xsh_polynomial_shift(*variance, 2, mean_y);
1751
1752 /* Maybe here add a consistency check that the variance polynomial is
1753 positive at all input points */
1754 }
1755
1756 /* If requested, compute mean squared error */
1757 if (mse != NULL || red_chisq != NULL) {
1758 int i;
1759
1760 if (mse != NULL)
1761 *mse = 0.00;
1762 if (red_chisq != NULL)
1763 *red_chisq = 0.00;
1764 for (i = 0; i < np; i++) {
1765 double regress = xsh_polynomial_evaluate_2d(result, xy_pos_data_x[i],
1766 xy_pos_data_y[i]);
1767 double residual = values_data[i] - regress;
1768 /* Subtract from the true value, square, accumulate */
1769 if (mse != NULL) {
1770 *mse += residual * residual;
1771 }
1772 if (red_chisq != NULL) {
1773 *red_chisq += (residual / sigmas_data[i]) * (residual / sigmas_data[i]) ;
1774 }
1775 }
1776 /* Get average */
1777 if (mse != NULL)
1778 *mse /= (double) np;
1779
1780 if (red_chisq != NULL) {
1781 passure( np > nc, "%d %d", np, nc);
1782 /* Was already checked */
1783 *red_chisq /= (double) (np - nc);
1784 }
1785 }
1786
1787 cleanup: return result;
1788}
1789
1790
static int dimension
static int degree
const char * xsh_tostring_cpl_type(cpl_type t)
Convert a CPL type to a string.
Definition: xsh_dump.c:359
#define assure(CONDITION, ERROR_CODE,...)
Definition: xsh_error.h:54
#define assure_nomsg(BOOL, CODE)
Definition: xsh_error.h:58
#define assure_mem(PTR)
Definition: xsh_error.h:79
#define check_msg(COMMAND,...)
Definition: xsh_error.h:62
#define passure(CONDITION,...)
Definition: xsh_error.h:82
int * x
#define xsh_msg_debug(...)
Print a debug message.
Definition: xsh_msg.h:99
polynomial * xsh_polynomial_new_zero(int dim)
Create a zero polynomial.
polynomial * xsh_polynomial_fit_2d(const cpl_bivector *xy_pos, const cpl_vector *values, const cpl_vector *sigmas, int poly_deg1, int poly_deg2, double *mse, double *red_chisq, polynomial **variance)
Fit a 2d surface with a polynomial in x and y.
cpl_matrix * xsh_matrix_solve_normal(const cpl_matrix *coeff, const cpl_matrix *rhs)
void xsh_polynomial_delete(polynomial **p)
Delete a polynomial.
double xsh_polynomial_evaluate_1d(const polynomial *p, double x)
Evaluate a 1d polynomial.
double xsh_polynomial_solve_1d(const polynomial *p, double value, double guess, int multiplicity)
Solve p(x) = value.
static cpl_error_code derivative_cpl_polynomial(cpl_polynomial *p, int varno)
Calculate the partial derivative of a CPL-polynomial.
polynomial * xsh_polynomial_collapse(const polynomial *p, int varno, double value)
Collapse a polynomial by fixing one variable to a constant.
void xsh_polynomial_delete_const(const polynomial **p)
Delete a const polynomial.
double xsh_polynomial_get_coeff_1d(const polynomial *p, int degree)
Get a coefficient of a 1D polynomial.
double xsh_polynomial_solve_2d(const polynomial *p, double value, double guess, int multiplicity, int varno, double x_value)
Solve p(x1, x2) = value.
polynomial * xsh_polynomial_fit_1d(const cpl_vector *x_pos, const cpl_vector *values, const cpl_vector *sigmas, int poly_deg, double *mse)
Fit a 1d function with a polynomial.
cpl_error_code xsh_polynomial_rescale(polynomial *p, int varno, double scale)
Rescale a polynomial.
double xsh_polynomial_derivative_1d(const polynomial *p, double x)
Evaluate the derivative of a 1d polynomial.
polynomial * xsh_polynomial_new(const cpl_polynomial *pol)
Create a polynomial.
polynomial * xsh_polynomial_add_2d(const polynomial *p1, const polynomial *p2)
Add two polynomials.
int xsh_polynomial_get_dimension(const polynomial *p)
Get the dimension of a polynomial.
int xsh_polynomial_get_degree(const polynomial *p)
Get degree.
polynomial * xsh_polynomial_convert_from_table(cpl_table *t)
Convert a table to a polynomial.
double xsh_polynomial_get_coeff_2d(const polynomial *p, int degree1, int degree2)
Get a coefficient of a 2D polynomial.
double xsh_polynomial_derivative_2d(const polynomial *p, double x1, double x2, int varno)
Evaluate the partial derivative of a 2d polynomial.
void xsh_polynomial_dump(const polynomial *p, FILE *stream)
Print a polynomial.
cpl_error_code xsh_polynomial_derivative(polynomial *p, int varno)
Calculate the partial derivative of a polynomial.
cpl_error_code xsh_polynomial_shift(polynomial *p, int varno, double shift)
Shift a polynomial.
cpl_matrix * xsh_matrix_product_normal_create(const cpl_matrix *self)
cpl_table * xsh_polynomial_convert_to_table(const polynomial *p)
Convert a polynomial to a table.
double xsh_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
Evaluate a 2d polynomial.
polynomial * xsh_polynomial_duplicate(const polynomial *p)
Copy a polynomial.
int xsh_max_int(int x, int y)
Maximum of two numbers.
Definition: xsh_utils.c:4415
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_matrix(cpl_matrix **m)
Deallocate a matrix and set the pointer to NULL.
Definition: xsh_utils.c:2209
double xsh_pow_int(double x, int y)
Computes x^y.
Definition: xsh_utils.c:4463
void xsh_free(const void *mem)
Deallocate memory.
Definition: xsh_utils.c:2102
cpl_polynomial * pol
cpl_vector * vec
int m
Definition: xsh_detmon_lg.c:91
int n
Definition: xsh_detmon_lg.c:92
DOUBLE mat[4][4]
#define COLUMN_ORDER1
#define COLUMN_COEFF
#define COLUMN_ORDER2