UVES Pipeline Reference Manual  5.4.6
uves_utils_polynomial.c
1 /* *
2  * This file is part of the ESO UVES 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: 2012-01-12 16:44:43 $
23  * $Revision: 1.68 $
24  * $Name: not supported by cvs2svn $
25  * $Log: not supported by cvs2svn $
26  * Revision 1.67 2011/12/08 14:03:32 amodigli
27  * Fix warnings with CPL6
28  *
29  * Revision 1.66 2010/09/24 09:32:08 amodigli
30  * put back QFITS dependency to fix problem spot by NRI on FIBER mode (with MIDAS calibs) data
31  *
32  * Revision 1.64 2007/09/11 17:08:49 amodigli
33  * mooved uves_polynomial_convert_from_plist_midas to uves_dfs
34  *
35  * Revision 1.63 2007/08/21 13:08:26 jmlarsen
36  * Removed irplib_access module, largely deprecated by CPL-4
37  *
38  * Revision 1.62 2007/06/20 15:34:50 jmlarsen
39  * Changed indentation
40  *
41  * Revision 1.61 2007/06/20 08:30:00 amodigli
42  * added index parameter to support FIBER mode lintab in uves_polynomial_convert_from_plist_midas
43  *
44  * Revision 1.60 2007/06/06 08:17:33 amodigli
45  * replace tab with 4 spaces
46  *
47  * Revision 1.59 2007/05/03 15:23:08 jmlarsen
48  * Removed dead code
49  *
50  * Revision 1.58 2007/05/03 15:18:29 jmlarsen
51  * Added function to add polynomials
52  *
53  * Revision 1.57 2007/04/27 07:21:51 jmlarsen
54  * Polyfit: Changed error code from ILLEGAL_INPUT to SINGULAR_MATRIX
55  *
56  * Revision 1.56 2007/04/24 12:50:29 jmlarsen
57  * Replaced cpl_propertylist -> uves_propertylist which is much faster
58  *
59  * Revision 1.55 2007/03/23 08:01:55 jmlarsen
60  * Fixed mixed code and declarations
61  *
62  * Revision 1.54 2007/03/19 15:10:03 jmlarsen
63  * Optimization in 2d fitting: do not call pow too often
64  *
65  * Revision 1.53 2007/03/13 15:35:11 jmlarsen
66  * Made a few time optimizations
67  *
68  * Revision 1.52 2007/03/05 10:20:49 jmlarsen
69  * Added uves_polynomial_delete_const()
70  *
71  * Revision 1.51 2007/01/15 08:48:51 jmlarsen
72  * Shortened lines
73  *
74  * Revision 1.50 2006/11/24 09:36:49 jmlarsen
75  * Workaround for slow uves_propertylist_get_size
76  *
77  * Revision 1.49 2006/11/15 15:02:15 jmlarsen
78  * Implemented const safe workarounds for CPL functions
79  *
80  * Revision 1.47 2006/11/15 14:04:08 jmlarsen
81  * Removed non-const version of parameterlist_get_first/last/next which is
82  * already in CPL, added const-safe wrapper, unwrapper and deallocator functions
83  *
84  * Revision 1.46 2006/11/13 14:23:55 jmlarsen
85  * Removed workarounds for CPL const bugs
86  *
87  * Revision 1.45 2006/11/06 15:19:42 jmlarsen
88  * Removed unused include directives
89  *
90  * Revision 1.44 2006/09/08 14:06:29 jmlarsen
91  * Removed profiling code
92  *
93  * Revision 1.43 2006/09/06 14:46:21 jmlarsen
94  * Added missing newline in uves_polynomial_dump()
95  *
96  * Revision 1.42 2006/08/17 14:11:25 jmlarsen
97  * Use assure_mem macro to check for memory allocation failure
98  *
99  * Revision 1.41 2006/08/17 13:56:53 jmlarsen
100  * Reduced max line length
101  *
102  * Revision 1.40 2006/07/03 13:27:52 jmlarsen
103  * Moved failing assertion to where it should be
104  *
105  * Revision 1.39 2006/06/01 14:43:17 jmlarsen
106  * Added missing documentation
107  *
108  * Revision 1.38 2006/05/05 13:59:03 jmlarsen
109  * Support fitting zero-degree polynomial
110  *
111  * Revision 1.37 2006/04/24 09:28:29 jmlarsen
112  * Added function to create zero-polynomial
113  *
114  * Revision 1.36 2006/03/27 09:41:37 jmlarsen
115  * Added timing markers
116  *
117  * Revision 1.35 2006/03/09 10:52:25 jmlarsen
118  * Renamed pow->power
119  *
120  * Revision 1.34 2006/03/03 13:54:11 jmlarsen
121  * Changed syntax of check macro
122  *
123  * Revision 1.33 2005/12/19 16:17:56 jmlarsen
124  * Replaced bool -> int
125  *
126  */
127 
128 #ifdef HAVE_CONFIG_H
129 # include <config.h>
130 #endif
131 
132 /*----------------------------------------------------------------------------*/
152 /*----------------------------------------------------------------------------*/
153 
154 /*-----------------------------------------------------------------------------
155  Defines
156  -----------------------------------------------------------------------------*/
157 
158 /*
159  * When storing a 2d polynomial in a table,
160  * these column names are used
161  */
162 #define COLUMN_ORDER1 "Order1"
163 #define COLUMN_ORDER2 "Order2"
164 #define COLUMN_COEFF "Coeff"
165 
168 /*-----------------------------------------------------------------------------
169  Includes
170  -----------------------------------------------------------------------------*/
171 #include <uves_utils_polynomial.h>
172 
173 #include <uves_utils.h>
174 #include <uves_utils_wrappers.h>
175 #include <uves_dump.h>
176 #include <uves_msg.h>
177 #include <uves_error.h>
178 
179 #include <cpl.h>
180 
181 /*-----------------------------------------------------------------------------
182  Typedefs
183  -----------------------------------------------------------------------------*/
186 struct _polynomial
187 {
189  cpl_polynomial *pol;
190 
192  cpl_vector *vec;
193  double *vec_data;
194 
195  int dimension; /* for efficiency */
196 
198  double *shift;
199 
201  double *scale;
202 };
203 
204 /*-----------------------------------------------------------------------------
205  Implementation
206  -----------------------------------------------------------------------------*/
207 /*----------------------------------------------------------------------------*/
218 /*----------------------------------------------------------------------------*/
219 polynomial *
220 uves_polynomial_new(const cpl_polynomial *pol)
221 {
222  polynomial *p = NULL;
223  int i;
224 
225  /* Test input */
226  assure(pol != NULL, CPL_ERROR_ILLEGAL_INPUT, "Null polynomial");
227 
228  /* Allocate and initialize struct */
229  p = cpl_calloc(1, sizeof(polynomial)) ;
230  assure_mem( p );
231 
232  check( p->dimension = cpl_polynomial_get_dimension(pol), "Error reading dimension");
233 
234  /* Allocate vector */
235  p->vec = cpl_vector_new(p->dimension);
236  assure_mem( p->vec );
237  p->vec_data = cpl_vector_get_data(p->vec);
238 
239  /* Shifts are initialized to zero, scales to 1 */
240  p->shift = cpl_calloc(p->dimension + 1, sizeof(double));
241  assure_mem( p->shift );
242 
243  p->scale = cpl_malloc((p->dimension + 1) * sizeof(double));
244  assure_mem( p->scale );
245  for (i = 0; i <= p->dimension; i++)
246  p->scale[i] = 1.0;
247 
248  check( p->pol = cpl_polynomial_duplicate(pol), "Error copying polynomial");
249 
250  cleanup:
251  if (cpl_error_get_code() != CPL_ERROR_NONE)
253 
254  return p;
255 }
256 
257 /*----------------------------------------------------------------------------*/
265 /*----------------------------------------------------------------------------*/
266 polynomial *
268 {
269  polynomial *result = NULL;
270  cpl_polynomial *p = NULL;
271 
272  assure( dim >= 1, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dim);
273 
274  p = cpl_polynomial_new(dim);
275  assure_mem( p );
276 
277  result = uves_polynomial_new(p);
278  assure_mem( result );
279 
280  cleanup:
281  uves_free_polynomial(&p);
282 
283  return result;
284 }
285 
286 /*----------------------------------------------------------------------------*/
293 /*----------------------------------------------------------------------------*/
294 void
296 {
298 }
299 
300 /*----------------------------------------------------------------------------*/
307 /*----------------------------------------------------------------------------*/
308 void
310 {
311  if (*p == NULL) return;
312  cpl_polynomial_delete((*p)->pol);
313  cpl_vector_delete((*p)->vec);
314  cpl_free((*p)->shift);
315  cpl_free((*p)->scale);
316  uves_free(*p);
317  *p = NULL;
318  return;
319 }
320 /*----------------------------------------------------------------------------*/
326 /*----------------------------------------------------------------------------*/
327 int
329 {
330  int result = -1;
331  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
332 
333  result = cpl_polynomial_get_degree(p->pol);
334 
335  cleanup:
336  return result;
337 }
338 
339 /*----------------------------------------------------------------------------*/
345 /*----------------------------------------------------------------------------*/
346 polynomial *
348 {
349  polynomial *result = NULL;
350  int dimension;
351  int i;
352 
353  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
354  dimension = uves_polynomial_get_dimension(p);
355 
356  check( result = uves_polynomial_new(p->pol),
357  "Error allocating polynomial");
358 
359  for (i = 0; i <= dimension; i++)
360  {
361  result->shift[i] = p->shift[i];
362  result->scale[i] = p->scale[i];
363  }
364 
365  cleanup:
366  if (cpl_error_get_code() != CPL_ERROR_NONE)
367  {
368  uves_polynomial_delete(&result);
369  return NULL;
370  }
371 
372  return result;
373 }
374 
375 
376 /*----------------------------------------------------------------------------*/
387 /*----------------------------------------------------------------------------*/
388 cpl_table *
390 {
391  cpl_table *t = NULL; /* Result */
392  int degree;
393  int i, j, row;
394 
395  /* Check input */
396  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
397  assure( uves_polynomial_get_dimension(p) == 2,
398  CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2D");
399 
400  degree = cpl_polynomial_get_degree(p->pol);
401 
402  /* Allocate space for 3 shifts, 3 scale factors and all
403  coefficients */
404  t = cpl_table_new(3 + 3 + (degree + 1)*(degree + 2)/2);
405  cpl_table_new_column(t, COLUMN_ORDER1, CPL_TYPE_INT);
406  cpl_table_new_column(t, COLUMN_ORDER2, CPL_TYPE_INT);
407  cpl_table_new_column(t, COLUMN_COEFF , CPL_TYPE_DOUBLE);
408 
409 
410  row = 0;
411 
412  /* First write the shifts, write non-garbage to coeff columns (which are not used) */
413  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
414  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
415  cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[0]); row++;
416 
417  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
418  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
419  cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[1]); row++;
420 
421  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
422  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
423  cpl_table_set_double(t, COLUMN_COEFF , row, p->shift[2]); row++;
424 
425  /* Then the scale factors */
426  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
427  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
428  cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[0]); row++;
429 
430  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
431  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
432  cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[1]); row++;
433 
434  cpl_table_set_int (t, COLUMN_ORDER1, row, -1);
435  cpl_table_set_int (t, COLUMN_ORDER2, row, -1);
436  cpl_table_set_double(t, COLUMN_COEFF, row, p->scale[2]); row++;
437 
438  /* And then write the coefficients */
439  for (i = 0; i <= degree; i++){
440  for (j = 0; j+i <= degree; j++){
441  double coeff;
442  cpl_size power[2];
443  power[0] = i;
444  power[1] = j;
445 
446  coeff = cpl_polynomial_get_coeff(p->pol, power);
447  cpl_table_set_int (t, COLUMN_ORDER1, row, power[0]);
448  cpl_table_set_int (t, COLUMN_ORDER2, row, power[1]);
449  cpl_table_set_double(t, COLUMN_COEFF , row, coeff);
450 
451  row++;
452  }
453  }
454  cpl_table_set_column_unit(t, COLUMN_ORDER1 , " ");
455  cpl_table_set_column_unit(t, COLUMN_ORDER2 , " ");
456  cpl_table_set_column_unit(t, COLUMN_COEFF , " ");
457 
458  cleanup:
459  return t;
460 }
461 
462 /*----------------------------------------------------------------------------*/
471 /*----------------------------------------------------------------------------*/
472 polynomial *
474 {
475  polynomial *p = NULL; /* Result */
476  cpl_polynomial *pol = NULL;
477  cpl_type type;
478  int i;
479 
480  /* Only 2d supported */
481  check( pol = cpl_polynomial_new(2), "Error initializing polynomial");
482 
483  /* Check table format */
484  assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
485  assure(cpl_table_has_column(t, COLUMN_ORDER1), CPL_ERROR_ILLEGAL_INPUT,
486  "No '%s' column found in table", COLUMN_ORDER1);
487  assure(cpl_table_has_column(t, COLUMN_ORDER2), CPL_ERROR_ILLEGAL_INPUT,
488  "No '%s' column found in table", COLUMN_ORDER2);
489  assure(cpl_table_has_column(t, COLUMN_COEFF ), CPL_ERROR_ILLEGAL_INPUT,
490  "No '%s' column found in table", COLUMN_COEFF );
491 
492  type = cpl_table_get_column_type(t, COLUMN_ORDER1);
493  assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
494  "Column '%s' has type %s. Integer expected", COLUMN_ORDER1,
495  uves_tostring_cpl_type(type));
496 
497  type = cpl_table_get_column_type(t, COLUMN_ORDER2);
498  assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
499  "Column '%s' has type %s. Integer expected", COLUMN_ORDER2,
500  uves_tostring_cpl_type(type));
501 
502  type = cpl_table_get_column_type(t, COLUMN_COEFF);
503  assure(type == CPL_TYPE_DOUBLE, CPL_ERROR_INVALID_TYPE,
504  "Column '%s' has type %s. Double expected", COLUMN_COEFF ,
505  uves_tostring_cpl_type(type));
506 
507  assure(cpl_table_get_nrow(t) > 1 + 2 + 1 + 2, CPL_ERROR_ILLEGAL_INPUT,
508  "Table must contain at least one coefficient");
509 
510  /* Read the coefficients */
511  for(i = 3 + 3; i < cpl_table_get_nrow(t); i++) {
512  double coeff;
513  cpl_size power[2];
514 
515  check(( power[0] = cpl_table_get_int(t, COLUMN_ORDER1, i, NULL),
516  power[1] = cpl_table_get_int(t, COLUMN_ORDER2, i, NULL),
517  coeff = cpl_table_get_double(t, COLUMN_COEFF , i, NULL)),
518  "Error reading table row %d", i);
519 
520  uves_msg_debug("Pol.coeff.(%" CPL_SIZE_FORMAT ", %" CPL_SIZE_FORMAT ") = %e", power[0], power[1], coeff);
521 
522  check( cpl_polynomial_set_coeff(pol, power, coeff), "Error creating polynomial");
523  }
524  p = uves_polynomial_new(pol);
525 
526  /* Read shifts and rescaling */
527  uves_polynomial_rescale(p, 0, cpl_table_get_double( t, COLUMN_COEFF, 3, NULL));
528  uves_polynomial_rescale(p, 1, cpl_table_get_double( t, COLUMN_COEFF, 4, NULL));
529  uves_polynomial_rescale(p, 2, cpl_table_get_double( t, COLUMN_COEFF, 5, NULL));
530  uves_polynomial_shift (p, 0, cpl_table_get_double( t, COLUMN_COEFF, 0, NULL));
531  uves_polynomial_shift (p, 1, cpl_table_get_double( t, COLUMN_COEFF, 1, NULL));
532  uves_polynomial_shift (p, 2, cpl_table_get_double( t, COLUMN_COEFF, 2, NULL));
533 
534  cleanup:
535  uves_free_polynomial(&pol);
536  if (cpl_error_get_code() != CPL_ERROR_NONE)
538 
539  return p;
540 }
541 
542 
543 /*----------------------------------------------------------------------------*/
549 /*----------------------------------------------------------------------------*/
550 int
552 {
553  int dim = -1;
554  assure(p != NULL, CPL_ERROR_ILLEGAL_INPUT, "Null polynomial");
555 
556 /* slow check( dim = cpl_polynomial_get_dimension(p->pol), "Error reading dimension"); */
557  dim = p->dimension;
558 
559  cleanup:
560  return dim;
561 }
562 
563 /*----------------------------------------------------------------------------*/
571 /*----------------------------------------------------------------------------*/
572 void uves_polynomial_dump(const polynomial *p, FILE *stream)
573 {
574  if (p == NULL)
575  fprintf(stream, "Null polynomial\n");
576  else {
577  int i;
578  cpl_polynomial_dump(p->pol, stream);
579  fprintf(stream, "shift_y \t= %f \tscale_y \t= %f\n", p->shift[0], p->scale[0]);
580  for (i = 1; i <= uves_polynomial_get_dimension(p); i++)
581  {
582  fprintf(stream, "shift_x%d \t= %f \tscale_x%d \t= %f\n",
583  i, p->shift[i], i, p->scale[i]);
584  }
585  }
586  return;
587 }
588 
589 /*----------------------------------------------------------------------------*/
603 /*----------------------------------------------------------------------------*/
604 cpl_error_code
605 uves_polynomial_rescale(polynomial *p, int varno, double scale)
606 {
607  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
608  assure(0 <= varno && varno <= uves_polynomial_get_dimension(p),
609  CPL_ERROR_ILLEGAL_INPUT, "Illegal variable number: %d", varno);
610 
611  /* Rescaling an x variable by the factor S corresponds to:
612  * p'(x) := p(x/S) =
613  * cpl( (x/S - shiftx ) / scalex ) * scaley + shifty =
614  * cpl( (x - (S*shiftx)) / (S*scalex) ) * scaley + shifty */
615 
616  /* Rescaling the y variable by the factor S corresponds to:
617  * p'(x) := S*p(x) =
618  * S * ( cpl((x - shiftx)/scalex) * scaley + shifty ) =
619  * cpl((x - shiftx)/scalex) * (S*scaley) + (S*shifty)
620  *
621  * therefore the implementation is the same in the two cases. */
622 
623  p->shift[varno] *= scale;
624  p->scale[varno] *= scale;
625 
626  cleanup:
627  return cpl_error_get_code();
628 }
629 
630 /*----------------------------------------------------------------------------*/
644 /*----------------------------------------------------------------------------*/
645 cpl_error_code
646 uves_polynomial_shift(polynomial *p, int varno, double shift)
647 {
648  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
649  assure(0 <= varno && varno <= uves_polynomial_get_dimension(p),
650  CPL_ERROR_ILLEGAL_INPUT, "Illegal variable number: %d", varno);
651 
652  /* The implementation is similar for x and y variables because
653  * p(x-S) =
654  * cpl( (x-S - shiftx) / scalex ) * scaley + shifty =
655  * cpl( (x - (shiftx+S)) / scalex ) * scaley + shifty
656  * and
657  * p(x) + S =
658  * cpl( (x - shiftx)/scalex ) * scaley + shifty + S =
659  * cpl( (x - shiftx)/scalex ) * scaley + (shifty+S) */
660 
661  p->shift[varno] += shift;
662 
663  cleanup:
664  return cpl_error_get_code();
665 }
666 
667 /*----------------------------------------------------------------------------*/
676 /*----------------------------------------------------------------------------*/
677 double
679 {
680  double result = 0;
681 
682  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
683  assure(uves_polynomial_get_dimension(p) == 1,
684  CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 1d");
685 
686  check( result =
687  cpl_polynomial_eval_1d(p->pol, (x - p->shift[1])/p->scale[1], NULL)
688  * p->scale[0] + p->shift[0],
689  "Could not evaluate polynomial");
690 
691  cleanup:
692  return result;
693 }
694 
695 
696 /*----------------------------------------------------------------------------*/
706 /*----------------------------------------------------------------------------*/
707 
708 double
709 uves_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
710 {
711  double result = 0;
712 
713  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
714  assure(p->dimension == 2, CPL_ERROR_ILLEGAL_INPUT,
715  "Polynomial must be 2d. It's %dd", p->dimension);
716  {
717  double scale = p->scale[0];
718  double shift = p->shift[0];
719 
720  // cpl_vector_set(p->vec, 0, (x1 - p->shift[1]) / p->scale[1]);
721  // cpl_vector_set(p->vec, 1, (x2 - p->shift[2]) / p->scale[2]);
722  p->vec_data[0] = (x1 - p->shift[1]) / p->scale[1];
723  p->vec_data[1] = (x2 - p->shift[2]) / p->scale[2];
724 
725  result = cpl_polynomial_eval(p->pol, p->vec) * scale + shift;
726  }
727 
728  cleanup:
729  return result;
730 }
731 
732 /*----------------------------------------------------------------------------*/
745 /*----------------------------------------------------------------------------*/
746 double
747 uves_polynomial_solve_1d(const polynomial *p, double value, double guess, int multiplicity)
748 {
749  double result = 0;
750  cpl_size power[1];
751  double coeff0;
752 
753  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
754  assure(uves_polynomial_get_dimension(p) == 1, CPL_ERROR_ILLEGAL_INPUT,
755  "Polynomial must be 1d");
756 
757  /* Solving p(x) = value corresponds to solving
758  <=> cpl_p( (x-xshift)/xscale )*yscale + yshift = value
759  <=> cpl_p( (x-xshift)/xscale ) + (yshift - value)/yscale = 0
760 
761  So 1) find zero point for the polynomial cpl() + (yshift-value)/yscale
762  Then 2) shift and rescale the result
763  */
764 
765  power[0] = 0;
766  check(( coeff0 = cpl_polynomial_get_coeff(p->pol, power),
767  cpl_polynomial_set_coeff(p->pol, power, coeff0 + (p->shift[0] - value)/p->scale[0])),
768  "Error setting coefficient");
769 
770  check( cpl_polynomial_solve_1d(p->pol, (guess - p->shift[1]) / p->scale[1],
771  &result, multiplicity), "Could not find root");
772  /* Restore polynomial */
773  cpl_polynomial_set_coeff(p->pol, power, coeff0);
774 
775  /* Shift solution */
776  result = result * p->scale[1] + p->shift[1];
777 
778  cleanup:
779  return result;
780 }
781 
782 /*----------------------------------------------------------------------------*/
799 /*----------------------------------------------------------------------------*/
800 double
801 uves_polynomial_solve_2d(const polynomial *p, double value, double guess,
802  int multiplicity, int varno, double x_value)
803 {
804  double result = 0;
805  polynomial *pol_1d = NULL;
806 
807  assure( 1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
808  "Illegal variable number: %d", varno);
809 
810  check( pol_1d = uves_polynomial_collapse(p, varno, x_value),
811  "Could not collapse polynomial");
812 
813  check( result = uves_polynomial_solve_1d(pol_1d, value, guess, multiplicity),
814  "Could not find root");
815 
816  cleanup:
817  uves_polynomial_delete(&pol_1d);
818  return result;
819 }
820 
821 /*----------------------------------------------------------------------------*/
830 /*----------------------------------------------------------------------------*/
831 double
832 uves_polynomial_derivative_2d(const polynomial *p, double x1, double x2, int varno)
833 {
834  double result = 0;
835  cpl_size power[2];
836 
837  assure (1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
838  "Illegal variable number (%d)", varno);
839 
840  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
841  assure(uves_polynomial_get_dimension(p) == 2, CPL_ERROR_ILLEGAL_INPUT,
842  "Polynomial must be 2d. It's %dd", uves_polynomial_get_dimension(p));
843 
844  /* d/dx_i [ p(x) ] =
845  * d/dx_i [ cpl( (x - shiftx) / scalex ) * scaley + shifty ] =
846  * [ d(cpl)/dx_i ( (x - shiftx) / scalex ) * scaley ]
847  */
848 
849  /* Shift, scale (x1, x2) */
850  x1 = (x1 - p->shift[1])/p->scale[1];
851  x2 = (x2 - p->shift[2])/p->scale[2];
852 
853  /* Get derivative of cpl polynomial.
854  *
855  */
856  {
857  int degree = cpl_polynomial_get_degree(p->pol);
858  double yj = 1; /* y^j */
859  int i, j;
860 
861  result = 0;
862  for (j = 0, yj = 1;
863  j <= degree; j++,
864  yj *= (varno == 1) ? x2 : x1)
865  {
866  /* Proof by example (degree = 3): For each j account for these terms
867  * using Horner's rule:
868  *
869  * d/dx y^j * [ c_3j x^3 + c_2j x^2 + c_1j x^1 + c_0j ] =
870  *
871  * y^j * [ 3c_3j x^2 + 2c_2j x^1 + 1c_1j ] =
872  *
873  * y^j * [ ((3c_3j) x + 2c_2j) x + 1c_1j ]
874  */
875 
876  double sum = 0;
877  for (i = degree; i >= 1; i--)
878  {
879  double c_ij;
880 
881  power[0] = (varno == 1) ? i : j;
882  power[1] = (varno == 1) ? j : i;
883 
884  c_ij = cpl_polynomial_get_coeff(p->pol, power);
885 
886  sum += (i * c_ij);
887  if (i >= 2) sum *= (varno == 1) ? x1 : x2;
888  }
889 
890  /* Collect terms */
891  result += yj * sum;
892  }
893  }
894 
895  result *= p->scale[0];
896 
897 
898 /* Old code: This method (valid for varno = 2)
899  of getting the derivative of
900  the CPL polynomial is slow because of the call to
901  uves_polynomial_collapse()
902 
903  check( pol_1d = uves_polynomial_collapse(p, 1, x1);
904  dummy = cpl_polynomial_eval_1d(pol_1d->pol, (x2 - pol_1d->shift[1])/pol_1d->scale[1], &result),
905  "Error evaluating derivative");
906 */
907 
908  cleanup:
909  return result;
910 }
911 
912 /*----------------------------------------------------------------------------*/
919 /*----------------------------------------------------------------------------*/
920 double
922 {
923  double result = 0;
924  double dummy;
925 
926  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
927  assure(uves_polynomial_get_dimension(p) == 1,
928  CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 1d");
929 
930  check( dummy = cpl_polynomial_eval_1d(p->pol, (x - p->shift[1])/p->scale[1], &result),
931  "Error evaluating derivative");
932 
933  cleanup:
934  return result;
935 }
936 
937 /*----------------------------------------------------------------------------*/
944 /*----------------------------------------------------------------------------*/
945 polynomial *
947 {
948  polynomial *result = NULL;
949  cpl_polynomial *pol = NULL;
950 
951  assure(p1 != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
952  assure(p2 != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
953  assure(uves_polynomial_get_dimension(p1) == 2,
954  CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
955  assure(uves_polynomial_get_dimension(p2) == 2,
956  CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
957 
958  /* cpl_polynomial1((x - shift_x1)/scale_x1) * scale_y1 + shift_y1
959  +
960  cpl_polynomial2((x - shift_x2)/scale_x2) * scale_y2 + shift_y2
961  = ???
962  Not easy.
963 
964  Use brute force:
965  */
966 
967  {
968  int degree, i, j;
969 
970  degree = uves_max_int(uves_polynomial_get_degree(p1),
972 
973  pol = cpl_polynomial_new(2);
974  for (i = 0; i <= degree; i++)
975  for (j = 0; j <= degree; j++) {
976  double coeff1, coeff2;
977  cpl_size power[2];
978 
979  /* Simple: add coefficients of the same power */
980  coeff1 = uves_polynomial_get_coeff_2d(p1, i, j);
981  coeff2 = uves_polynomial_get_coeff_2d(p2, i, j);
982 
983  power[0] = i;
984  power[1] = j;
985  cpl_polynomial_set_coeff(pol, power, coeff1 + coeff2);
986  }
987  }
988 
989  result = uves_polynomial_new(pol);
990 
991  cleanup:
992  uves_free_polynomial(&pol);
993  return result;
994 }
995 
996 /*----------------------------------------------------------------------------*/
1009 /*----------------------------------------------------------------------------*/
1010 static cpl_error_code
1011 derivative_cpl_polynomial(cpl_polynomial *p, int varno)
1012 {
1013  int dimension, degree;
1014  int i, j;
1015  cpl_size power[2];
1016 
1017  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1018  dimension = cpl_polynomial_get_dimension(p);
1019  degree = cpl_polynomial_get_degree(p);
1020  assure( 1 <= dimension && dimension <= 2, CPL_ERROR_ILLEGAL_INPUT,
1021  "Illegal dimension: %d", dimension);
1022  assure( 1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1023  "Illegal variable number: %d", varno);
1024 
1025  if (dimension == 1)
1026  {
1027  /* a_i := (i+1) * a_(i+1) */
1028  for(i = 0; i <= degree; i++)
1029  {
1030  double coeff;
1031  power[0] = i+1;
1032  /* power[1] is ignored */
1033 
1034  coeff = cpl_polynomial_get_coeff(p, power);
1035 
1036  power[0] = i;
1037  cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1038  }
1039  }
1040 
1041  if (dimension == 2)
1042  {
1043  /* a_ij := (i+1) * a_{(i+1),j} */
1044  for(i = 0; i <= degree; i++)
1045  {
1046  for(j = 0; i + j <= degree; j++)
1047  {
1048  double coeff;
1049  power[varno - 1] = i+1; /* varno == 1: 0,1 */
1050  power[2 - varno] = j; /* varno == 2: 1,0 */
1051 
1052  coeff = cpl_polynomial_get_coeff(p, power);
1053 
1054  power[varno - 1] = i;
1055 
1056  cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1057  }
1058  }
1059  }
1060 
1061  cleanup:
1062  return cpl_error_get_code();
1063 }
1064 
1065 /*----------------------------------------------------------------------------*/
1075 /*----------------------------------------------------------------------------*/
1076 cpl_error_code
1078 {
1079  int dimension;
1080 
1081  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1082  check ( dimension = uves_polynomial_get_dimension(p), "Error reading dimension");
1083  assure( 1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1084  "Illegal variable number: %d", varno);
1085 
1086 
1087  /* d/dx_i [ cpl( (x - shiftx) / scalex ) * scaley + shifty ] =
1088  * sum_j d(cpl)/dx_j ( (x - shiftx) / scalex ) * scaley * dx_j/dx_i / scalex_j =
1089  * d(cpl)/dx_i ( (x - shiftx) / scalex ) * scaley/scalex_i,
1090  *
1091  * so transform : shifty -> 0
1092  * shiftx -> shiftx
1093  * scaley -> scaley/scalex_i
1094  * scalex -> scalex
1095  * cpl -> d(cpl)/dx_i
1096  */
1097 
1098  p->shift[0] = 0;
1099  p->scale[0] = p->scale[0] / p->scale[varno];
1100 
1101  check( derivative_cpl_polynomial(p->pol, varno),
1102  "Error calculating derivative of CPL-polynomial");
1103 
1104  cleanup:
1105  return cpl_error_get_code();
1106 }
1107 
1108 
1109 /*----------------------------------------------------------------------------*/
1118 /*----------------------------------------------------------------------------*/
1119 double
1120 uves_polynomial_get_coeff_2d(const polynomial *p, int degree1, int degree2)
1121 {
1122  polynomial *pp = NULL;
1123  int dimension;
1124  double result = 0;
1125  double factorial;
1126 
1127  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1128  check ( dimension = uves_polynomial_get_dimension(p), "Error reading dimension");
1129  assure(dimension == 2, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dimension);
1130  assure( 0 <= degree1, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree1);
1131  assure( 0 <= degree2, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree2);
1132 
1133  /* Calculate the coefficient as
1134  * d^N p / (dx1^degree1 dx2^degree2) / (degree1! * degree2!)
1135  * evaluated in (0,0)
1136  */
1137 
1138  pp = uves_polynomial_duplicate(p);
1139 
1140  factorial = 1;
1141  while(degree1 > 0)
1142  {
1143  check( uves_polynomial_derivative(pp, 1), "Error calculating derivative");
1144 
1145  factorial *= degree1;
1146  degree1 -= 1;
1147  }
1148 
1149  while(degree2 > 0)
1150  {
1151  check( uves_polynomial_derivative(pp, 2), "Error calculating derivative");
1152 
1153  factorial *= degree2;
1154  degree2 -= 1;
1155  }
1156 
1157  check( result = uves_polynomial_evaluate_2d(pp, 0, 0) / factorial,
1158  "Error evaluating polynomial");
1159 
1160  cleanup:
1162  return result;
1163 }
1164 /*----------------------------------------------------------------------------*/
1174 /*----------------------------------------------------------------------------*/
1175 double
1177 {
1178  polynomial *pp = NULL;
1179  int dimension;
1180  double result = 0;
1181  double factorial;
1182 
1183  assure( p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1184  check ( dimension = uves_polynomial_get_dimension(p), "Error reading dimension");
1185  assure(dimension == 1, CPL_ERROR_ILLEGAL_INPUT, "Illegal dimension: %d", dimension);
1186  assure( 0 <= degree, CPL_ERROR_ILLEGAL_INPUT, "Illegal degree: %d", degree);
1187 
1188  /* Calculate the coefficient as
1189  * d^degree p/dx^degree / (degree1! * degree2!)
1190  * evaluated in 0.
1191  */
1192 
1193  pp = uves_polynomial_duplicate(p);
1194 
1195  factorial = 1;
1196  while(degree > 0)
1197  {
1198  check( uves_polynomial_derivative(pp, 1), "Error calculating derivative");
1199 
1200  factorial *= degree;
1201  degree -= 1;
1202  }
1203 
1204  check( result = uves_polynomial_evaluate_1d(pp, 0) / factorial,
1205  "Error evaluating polynomial");
1206 
1207  cleanup:
1209  return result;
1210 }
1211 
1212 
1213 /*----------------------------------------------------------------------------*/
1229 /*----------------------------------------------------------------------------*/
1230 polynomial *
1231 uves_polynomial_collapse(const polynomial *p, int varno, double value)
1232 {
1233  polynomial *result = NULL;
1234  cpl_polynomial *pol = NULL;
1235  cpl_size *power = NULL;
1236 
1237  int i, j;
1238  int degree, dimension;
1239 
1240  assure(p != NULL, CPL_ERROR_NULL_INPUT, "Null polynomial");
1241  dimension = uves_polynomial_get_dimension(p);
1242  assure(dimension > 0, CPL_ERROR_ILLEGAL_INPUT,
1243  "Polynomial has non-positive dimension: %d", dimension);
1244  assure(dimension != 1, CPL_ERROR_ILLEGAL_OUTPUT,
1245  "Don't collapse a 1d polynomial. Evaluate it!");
1246 
1247  /* To generalize this function to work with dimensions higher than 2,
1248  also changes needs to be made below (use varno properly). For now,
1249  support only 2d. */
1250  assure(dimension == 2, CPL_ERROR_ILLEGAL_INPUT, "Polynomial must be 2d");
1251 
1252  assure(1 <= varno && varno <= dimension, CPL_ERROR_ILLEGAL_INPUT,
1253  "Wrong variable number");
1254  value = (value - p->shift[varno]) / p->scale[varno];
1255 
1256  /* Compute new coefficients */
1257  degree = cpl_polynomial_get_degree(p->pol);
1258  pol = cpl_polynomial_new(dimension - 1);
1259  power = cpl_malloc(sizeof(cpl_size) * dimension);
1260  assure_mem( power );
1261  for (i = 0; i <= degree; i++)
1262  {
1263  double coeff;
1264 
1265  power[2-varno] = i; /* map 2->0 and 1->1 */
1266 
1267  /* Collect all terms with x^i (using Horner's rule) */
1268  coeff = 0;
1269  for (j = degree - i; j >= 0; j--)
1270  {
1271  power[varno-1] = j; /* map 2->1 and 1->0 */
1272  coeff += cpl_polynomial_get_coeff(p->pol, power);
1273  if (j > 0) coeff *= value;
1274  }
1275  /* Write coefficient in 1d polynomial */
1276  power[0] = i;
1277  cpl_polynomial_set_coeff(pol, power, coeff);
1278  }
1279 
1280  /* Wrap the polynomial */
1281  result = uves_polynomial_new(pol);
1282 
1283  /* Copy the shifts and scales, skip variable number varno */
1284  j = 0;
1285  for(i = 0; i <= dimension - 1; i++)
1286  {
1287  if (i == varno)
1288  {
1289  /* Don't copy */
1290  j += 2;
1291  /* For the remainder of this for loop, j = i+1 */
1292  }
1293  else
1294  {
1295  result->shift[i] = p->shift[j];
1296  result->scale[i] = p->scale[j];
1297  j += 1;
1298  }
1299  }
1300 
1301  assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
1302  "Error collapsing polynomial");
1303 
1304  cleanup:
1305  cpl_free(power); power = NULL;
1306  uves_free_polynomial(&pol);
1307  if (cpl_error_get_code() != CPL_ERROR_NONE)
1308  {
1309  uves_polynomial_delete(&result);
1310  }
1311  return result;
1312 }
1313 
1314 
1315 
1316 /*----------------------------------------------------------------------------*/
1336 /*----------------------------------------------------------------------------*/
1338  const cpl_vector * x_pos,
1339  const cpl_vector * values,
1340  const cpl_vector * sigmas,
1341  int poly_deg,
1342  double * mse)
1343 {
1344  int nc ;
1345  int np ;
1346  cpl_matrix * ma = NULL;
1347  cpl_matrix * mb = NULL;
1348  cpl_matrix * mx = NULL;
1349  const double * x_pos_data ;
1350  const double * values_data ;
1351  const double * sigmas_data = NULL;
1352  double mean_x, mean_z;
1353  polynomial * result = NULL;
1354  cpl_polynomial * out ;
1355  cpl_vector * x_val = NULL;
1356  int i, j ;
1357 
1358  /* Check entries */
1359  assure_nomsg( x_pos != NULL && values != NULL, CPL_ERROR_NULL_INPUT);
1360  assure( poly_deg >= 0, CPL_ERROR_ILLEGAL_INPUT,
1361  "Polynomial degree is %d. Must be non-negative", poly_deg);
1362  np = cpl_vector_get_size(x_pos) ;
1363 
1364  nc = 1 + poly_deg ;
1365  assure( np >= nc, CPL_ERROR_ILLEGAL_INPUT,
1366  "Not enough points (%d) to fit %d-order polynomial. %d point(s) needed",
1367  np, poly_deg, nc);
1368 
1369  /* Fill up look-up table for coefficients to compute */
1370  /* Initialize matrices */
1371  /* ma contains the polynomial terms for each input point. */
1372  /* mb contains the values */
1373  ma = cpl_matrix_new(np, nc) ;
1374  mb = cpl_matrix_new(np, 1) ;
1375 
1376  /* Get mean values */
1377  mean_x = cpl_vector_get_mean(x_pos);
1378  mean_z = cpl_vector_get_mean(values);
1379 
1380  /* Fill up matrices, shift */
1381  x_pos_data = cpl_vector_get_data_const(x_pos) ;
1382  values_data = cpl_vector_get_data_const(values) ;
1383  if (sigmas != NULL)
1384  {
1385  sigmas_data = cpl_vector_get_data_const(sigmas) ;
1386  }
1387 
1388  if (sigmas != NULL)
1389  {
1390  for (i=0 ; i<np ; i++)
1391  {
1392  /* Catch division by zero */
1393  if (sigmas_data[i] == 0)
1394  {
1395  uves_free_matrix(&ma) ;
1396  uves_free_matrix(&mb) ;
1397  assure(false, CPL_ERROR_DIVISION_BY_ZERO,
1398  "Sigmas must be non-zero");
1399  }
1400  for (j=0 ; j<nc ; j++)
1401  {
1402  cpl_matrix_set(ma, i, j,
1403  uves_pow_int(x_pos_data[i] - mean_x, j) /
1404  sigmas_data[i]) ;
1405  }
1406  /* mb contains surface values (z-axis) */
1407  cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / sigmas_data[i]);
1408  }
1409  }
1410  else /* Use sigma = 1 */
1411  {
1412  for (i=0 ; i<np ; i++)
1413  {
1414  for (j=0 ; j<nc ; j++)
1415  {
1416  cpl_matrix_set(ma, i, j,
1417  uves_pow_int(x_pos_data[i] - mean_x, j) / 1);
1418  }
1419  /* mb contains surface values (z-values) */
1420  cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / 1) ;
1421  }
1422  }
1423 
1424  /* Solve XA=B by a least-square solution (aka pseudo-inverse). */
1425  check( mx = cpl_matrix_solve_normal(ma, mb),
1426  "Could not invert matrix");
1427  uves_free_matrix(&ma);
1428  uves_free_matrix(&mb);
1429 
1430  /* Store coefficients for output */
1431  out = cpl_polynomial_new(1) ;
1432  cpl_size deg=0;
1433  for (deg=0 ; deg<nc ; deg++) {
1434  cpl_polynomial_set_coeff(out, &deg, cpl_matrix_get(mx, deg, 0)) ;
1435  }
1436  uves_free_matrix(&mx);
1437 
1438  /* If requested, compute mean squared error */
1439  if (mse != NULL) {
1440  *mse = 0.00 ;
1441  x_val = cpl_vector_new(1) ;
1442  for (i=0 ; i<np ; i++)
1443  {
1444  double residual;
1445  cpl_vector_set(x_val, 0, x_pos_data[i] - mean_x) ;
1446  /* Subtract from the true value, square, accumulate */
1447  residual = (values_data[i] - mean_z) - cpl_polynomial_eval(out, x_val);
1448  *mse += residual*residual;
1449  }
1450  uves_free_vector(&x_val) ;
1451  /* Average the error term */
1452  *mse /= (double)np ;
1453  }
1454 
1455  /* Create and shift result */
1456  result = uves_polynomial_new(out);
1457  uves_free_polynomial(&out);
1458 
1459  uves_polynomial_shift(result, 0, mean_z);
1460  uves_polynomial_shift(result, 1, mean_x);
1461 
1462  cleanup:
1463  uves_free_vector(&x_val);
1464  uves_free_matrix(&ma);
1465  uves_free_matrix(&mb);
1466  uves_free_matrix(&mx);
1467  return result;
1468 }
1469 
1470 
1471 /*----------------------------------------------------------------------------*/
1515 /*----------------------------------------------------------------------------*/
1516 polynomial *
1518  const cpl_bivector * xy_pos,
1519  const cpl_vector * values,
1520  const cpl_vector * sigmas,
1521  int poly_deg1,
1522  int poly_deg2,
1523  double * mse,
1524  double * red_chisq,
1525  polynomial ** variance)
1526 {
1527  int nc ;
1528  int degx, degy ;
1529  int * degx_tab ;
1530  int * degy_tab ;
1531  int np ;
1532  cpl_matrix * ma ;
1533  cpl_matrix * mb ;
1534  cpl_matrix * mx ;
1535  cpl_matrix * mat;
1536  cpl_matrix * mat_ma;
1537  cpl_matrix * cov = NULL;
1538  const double * xy_pos_data_x ;
1539  const double * xy_pos_data_y ;
1540  const double * values_data ;
1541  const double * sigmas_data = NULL;
1542  const cpl_vector* xy_pos_x;
1543  const cpl_vector* xy_pos_y;
1544  double mean_x, mean_y, mean_z;
1545  cpl_polynomial * out ;
1546  cpl_polynomial * variance_cpl ;
1547  polynomial * result = NULL;
1548  cpl_size * powers ;
1549 
1550  /* Check entries */
1551  assure(xy_pos && values, CPL_ERROR_NULL_INPUT, "Null input");
1552  assure(poly_deg1 >= 0, CPL_ERROR_ILLEGAL_INPUT, "Polynomial degree1 is %d", poly_deg1);
1553  assure(poly_deg2 >= 0, CPL_ERROR_ILLEGAL_INPUT, "Polynomial degree2 is %d", poly_deg2);
1554  np = cpl_bivector_get_size(xy_pos) ;
1555 
1556  /* Can't calculate variance and chi_sq without sigmas */
1557  assure( (variance == NULL && red_chisq == NULL) || sigmas != NULL,
1558  CPL_ERROR_ILLEGAL_INPUT,
1559  "Cannot calculate variance or chi_sq without knowing");
1560 
1561  /* Fill up look-up table for coefficients to compute */
1562  nc = (1 + poly_deg1)*(1 + poly_deg2) ; /* rectangular matrix */
1563 
1564  assure(np >= nc, CPL_ERROR_SINGULAR_MATRIX, "%d coefficients. Only %d points", nc, np);
1565  /* The error code here is set to SINGULAR_MATRIX, in order to allow the caller
1566  to detect when too many coefficients are fitted to too few points */
1567 
1568  /* Need an extra point to calculate reduced chi^2 */
1569  assure(red_chisq == NULL || np > nc, CPL_ERROR_ILLEGAL_INPUT,
1570  "%d coefficients. %d points. Cannot calculate chi square", nc, np);
1571 
1572  degx_tab = cpl_malloc(nc * sizeof(int)) ;
1573  assure_mem( degx_tab );
1574 
1575  degy_tab = cpl_malloc(nc * sizeof(int)) ;
1576  if (degy_tab == NULL) {
1577  cpl_free(degx_tab);
1578  assure_mem( false );
1579  }
1580 
1581  {
1582  int i=0 ;
1583  for (degy=0 ; degy<=poly_deg2 ; degy++) { /* rectangular matrix */
1584  for (degx=0 ; degx<=poly_deg1 ; degx++) {
1585  degx_tab[i] = degx ;
1586  degy_tab[i] = degy ;
1587  i++ ;
1588  }
1589  }
1590  }
1591 
1592  /* Initialize matrices */
1593  /* ma contains the polynomial terms in the order described */
1594  /* above in each column, for each input point. */
1595  /* mb contains the values */
1596  ma = cpl_matrix_new(np, nc) ;
1597  mb = cpl_matrix_new(np, 1) ;
1598 
1599  /* Get the mean of each variable */
1600  xy_pos_x = cpl_bivector_get_x_const(xy_pos);
1601  xy_pos_y = cpl_bivector_get_y_const(xy_pos);
1602 
1603  mean_x = cpl_vector_get_mean(xy_pos_x);
1604  mean_y = cpl_vector_get_mean(xy_pos_y);
1605  mean_z = cpl_vector_get_mean(values);
1606 
1607  /* Fill up matrices. At the same time shift the data
1608  so that it is centered around zero */
1609  xy_pos_data_x = cpl_vector_get_data_const(xy_pos_x) ;
1610  xy_pos_data_y = cpl_vector_get_data_const(xy_pos_y) ;
1611  values_data = cpl_vector_get_data_const(values) ;
1612  if (sigmas != NULL)
1613  {
1614  sigmas_data = cpl_vector_get_data_const(sigmas) ;
1615  }
1616 
1617  if (sigmas != NULL)
1618  {
1619  int i;
1620  for (i=0 ; i<np ; i++) {
1621  double *ma_data = cpl_matrix_get_data(ma);
1622  double *mb_data = cpl_matrix_get_data(mb);
1623 
1624  int j = 0;
1625  double valy = 1;
1626 
1627  /* Catch division by zero */
1628  if (sigmas_data[i] == 0)
1629  {
1630  uves_free_matrix(&ma) ;
1631  uves_free_matrix(&mb) ;
1632  cpl_free(degx_tab) ;
1633  cpl_free(degy_tab) ;
1634  assure(false, CPL_ERROR_DIVISION_BY_ZERO,
1635  "Sigmas must be non-zero. sigma[%d] is %f", i, sigmas_data[i]);
1636  }
1637 
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 / sigmas_data[i];
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 
1650  mb_data[0 + i*1] = (values_data[i] - mean_z) / sigmas_data[i];
1651  }
1652  }
1653  else /* Use sigma = 1 */
1654  {
1655  int i;
1656  for (i=0 ; i<np ; i++) {
1657  double *ma_data = cpl_matrix_get_data(ma);
1658  double *mb_data = cpl_matrix_get_data(mb);
1659 
1660  double valy = 1;
1661  int j = 0;
1662  for (degy=0 ; degy<=poly_deg2 ; degy++) {
1663  double valx = 1;
1664  for (degx=0 ; degx<=poly_deg1 ; degx++) {
1665  ma_data[j + i*nc] = valx * valy / 1;
1666  valx *= (xy_pos_data_x[i] - mean_x);
1667  j++;
1668  }
1669  valy *= (xy_pos_data_y[i] - mean_y);
1670  }
1671 
1672  /* mb contains surface values (z-axis) */
1673 // cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / 1) ;
1674  mb_data[0 + i*1] = (values_data[i] - mean_z) / 1;
1675  }
1676  }
1677 
1678  /* If variance polynomial is requested,
1679  compute covariance matrix = (A^T * A)^-1 */
1680  if (variance != NULL)
1681  {
1682  mat = cpl_matrix_transpose_create(ma);
1683  if (mat != NULL)
1684  {
1685  mat_ma = cpl_matrix_product_create(mat, ma);
1686  if (mat_ma != NULL)
1687  {
1688  cov = cpl_matrix_invert_create(mat_ma);
1689  /* Here, one might do a (paranoia) check that the covariance
1690  matrix is symmetrical and has positive eigenvalues (so that
1691  the returned variance polynomial is guaranteed to be positive) */
1692 
1693  variance_cpl = cpl_polynomial_new(2);
1694  }
1695  }
1696  uves_free_matrix(&mat);
1697  uves_free_matrix(&mat_ma);
1698  }
1699 
1700  /* Solve XA=B by a least-square solution (aka pseudo-inverse). */
1701  mx = cpl_matrix_solve_normal(ma, mb) ;
1702 
1703  uves_free_matrix(&ma) ;
1704  uves_free_matrix(&mb) ;
1705  if (mx == NULL) {
1706  cpl_free(degx_tab) ;
1707  cpl_free(degy_tab) ;
1708  uves_free_matrix(&cov) ;
1709  assure(false, CPL_ERROR_ILLEGAL_OUTPUT, "Matrix inversion failed") ;
1710  }
1711 
1712  /* Store coefficients for output */
1713  out = cpl_polynomial_new(2) ;
1714  powers = cpl_malloc(2 * sizeof(cpl_size)) ;
1715  if (powers == NULL) {
1716  cpl_free(degx_tab) ;
1717  cpl_free(degy_tab) ;
1718  uves_free_matrix(&mx) ;
1719  uves_free_matrix(&cov) ;
1720  uves_free_polynomial(&out) ;
1721  assure_mem( false );
1722  }
1723 
1724  {
1725  int i;
1726  for (i = 0 ; i < nc ; i++)
1727  {
1728  powers[0] = degx_tab[i] ;
1729  powers[1] = degy_tab[i] ;
1730  cpl_polynomial_set_coeff(out, powers, cpl_matrix_get(mx, i, 0)) ;
1731 
1732  /* Create variance polynomial (if requested) */
1733  if (variance != NULL && /* Requested? */
1734  cov != NULL && variance_cpl != NULL /* covariance computation succeeded? */
1735  )
1736  {
1737  int j;
1738  for (j = 0; j < nc; j++)
1739  {
1740  double coeff;
1741  /* Add cov_ij to the proper coeff:
1742  cov_ij * dp/d(ai) * dp/d(aj) =
1743  cov_ij * (x^degx[i] * y^degy[i]) * (x^degx[i] * y^degy[i]) =
1744  cov_ij * x^(degx[i]+degx[j]) * y^(degy[i] + degy[j]),
1745 
1746  i.e. add cov_ij to coeff (degx[i]+degx[j], degy[i]+degy[j]) */
1747  powers[0] = degx_tab[i] + degx_tab[j] ;
1748  powers[1] = degy_tab[i] + degy_tab[j] ;
1749 
1750  coeff = cpl_polynomial_get_coeff(variance_cpl, powers);
1751  cpl_polynomial_set_coeff(variance_cpl, powers,
1752  coeff + cpl_matrix_get(cov, i, j)) ;
1753  }
1754  }
1755  }
1756  }
1757 
1758  cpl_free(powers) ;
1759  cpl_free(degx_tab) ;
1760  cpl_free(degy_tab) ;
1761  uves_free_matrix(&cov) ;
1762  uves_free_matrix(&mx) ;
1763 
1764  /* Create and shift result */
1765  result = uves_polynomial_new(out);
1766  uves_free_polynomial(&out);
1767  uves_polynomial_shift(result, 0, mean_z);
1768  uves_polynomial_shift(result, 1, mean_x);
1769  uves_polynomial_shift(result, 2, mean_y);
1770 
1771  /* Wrap up variance polynomial */
1772  if (variance != NULL)
1773  {
1774  *variance = uves_polynomial_new(variance_cpl);
1775  uves_free_polynomial(&variance_cpl);
1776  /* The variance of the fit does not change
1777  when a constant is added to the a_00
1778  coefficient of the polynomial, so don't:
1779  uves_polynomial_shift(*variance, 0, mean_z); */
1780  uves_polynomial_shift(*variance, 1, mean_x);
1781  uves_polynomial_shift(*variance, 2, mean_y);
1782 
1783  /* Maybe here add a consistency check that the variance polynomial is
1784  positive at all input points */
1785  }
1786 
1787  /* If requested, compute mean squared error */
1788  if (mse != NULL || red_chisq != NULL)
1789  {
1790  int i;
1791 
1792  if (mse != NULL) *mse = 0.00 ;
1793  if (red_chisq != NULL) *red_chisq = 0.00 ;
1794  for (i = 0 ; i < np ; i++)
1795  {
1796  double regress = uves_polynomial_evaluate_2d(result,
1797  xy_pos_data_x[i],
1798  xy_pos_data_y[i]);
1799  /* Subtract from the true value, square, accumulate */
1800  if (mse != NULL)
1801  {
1802  double residual = values_data[i] - regress;
1803  *mse += residual*residual;
1804  }
1805  if (red_chisq != NULL)
1806  {
1807  *red_chisq += uves_pow_int((values_data[i] - regress) /
1808  sigmas_data[i], 2);
1809  }
1810  }
1811  /* Get average */
1812  if (mse != NULL) *mse /= (double) np ;
1813 
1814  if (red_chisq != NULL)
1815  {
1816  passure( np > nc, "%d %d", np, nc); /* Was already checked */
1817  *red_chisq /= (double) (np - nc) ;
1818  }
1819  }
1820 
1821  cleanup:
1822  return result ;
1823 }
1824 
1825 
int uves_polynomial_get_dimension(const polynomial *p)
Get the dimension of a polynomial.
polynomial * uves_polynomial_collapse(const polynomial *p, int varno, double value)
Collapse a polynomial by fixing one variable to a constant.
polynomial * uves_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.
double uves_polynomial_solve_2d(const polynomial *p, double value, double guess, int multiplicity, int varno, double x_value)
Solve p(x1, x2) = value.
double uves_polynomial_derivative_1d(const polynomial *p, double x)
Evaluate the derivative of a 1d polynomial.
void uves_polynomial_delete(polynomial **p)
Delete a polynomial.
cpl_error_code uves_polynomial_derivative(polynomial *p, int varno)
Calculate the partial derivative of a polynomial.
polynomial * uves_polynomial_add_2d(const polynomial *p1, const polynomial *p2)
Add two polynomials.
polynomial * uves_polynomial_convert_from_table(cpl_table *t)
Convert a table to a polynomial.
static cpl_error_code derivative_cpl_polynomial(cpl_polynomial *p, int varno)
Calculate the partial derivative of a CPL-polynomial.
double uves_pow_int(double x, int y)
Calculate x to the y'th.
Definition: uves_utils.c:1593
#define passure(BOOL,...)
Definition: uves_error.h:207
cpl_table * uves_polynomial_convert_to_table(const polynomial *p)
Convert a polynomial to a table.
cpl_polynomial * pol
double uves_polynomial_get_coeff_2d(const polynomial *p, int degree1, int degree2)
Get a coefficient of a 2D polynomial.
double uves_polynomial_derivative_2d(const polynomial *p, double x1, double x2, int varno)
Evaluate the partial derivative of a 2d polynomial.
int uves_polynomial_get_degree(const polynomial *p)
Get degree.
polynomial * uves_polynomial_duplicate(const polynomial *p)
Copy a polynomial.
void uves_polynomial_delete_const(const polynomial **p)
Delete a const polynomial.
cpl_error_code uves_polynomial_rescale(polynomial *p, int varno, double scale)
Rescale a polynomial.
#define assure_mem(PTR)
Definition: uves_error.h:181
polynomial * uves_polynomial_new(const cpl_polynomial *pol)
Create a polynomial.
double uves_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
Evaluate a 2d polynomial.
double uves_polynomial_evaluate_1d(const polynomial *p, double x)
Evaluate a 1d polynomial.
double uves_polynomial_get_coeff_1d(const polynomial *p, int degree)
Get a coefficient of a 1D polynomial.
void uves_polynomial_dump(const polynomial *p, FILE *stream)
Print a polynomial.
const char * uves_tostring_cpl_type(cpl_type t)
Convert a CPL type to a string.
Definition: uves_dump.c:378
cpl_error_code uves_polynomial_shift(polynomial *p, int varno, double shift)
Shift a polynomial.
#define uves_msg_debug(...)
Print a debug message.
Definition: uves_msg.h:97
#define assure_nomsg(BOOL, CODE)
Definition: uves_error.h:177
#define check(CMD,...)
Definition: uves_error.h:198
polynomial * uves_polynomial_new_zero(int dim)
Create a zero polynomial.
double uves_polynomial_solve_1d(const polynomial *p, double value, double guess, int multiplicity)
Solve p(x) = value.
polynomial * uves_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.