CR2RE Pipeline Reference Manual 1.6.10
hdrl_prototyping.c
1/*
2 * This file is part of the HDRL
3 * Copyright (C) 2014 European Southern Observatory
4 *
5 * This program 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, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20
21#ifdef HAVE_CONFIG_H
22#include <config.h>
23#endif
24
25/*-----------------------------------------------------------------------------
26 Includes
27 -----------------------------------------------------------------------------*/
28
29#include "hdrl_prototyping.h"
30
31#include <stdio.h>
32#include <stdlib.h>
33#include <math.h>
34#include <errno.h>
35#include <complex.h>
36/*----------------------------------------------------------------------------*/
44/*----------------------------------------------------------------------------*/
45
46/*-----------------------------------------------------------------------------
47 Defines
48 -----------------------------------------------------------------------------*/
49/*-----------------------------------------------------------------------------
50 Functions
51 -----------------------------------------------------------------------------*/
52
56/*-----------------------------------------------------------------------------
57 Static
58 -----------------------------------------------------------------------------*/
59
60
61static cpl_image *hdrl_mirror_edges(cpl_image *image, int dx, int dy);
62static cpl_image *hdrl_gen_lowpass(int xs, int ys, double sigma_x,
63 double sigma_y);
64
65
66
77cpl_image * hdrl_get_spatial_freq(cpl_image *ima, double gausfilt, int mirrorx,
78 int mirrory){
79
80 int xsize=0, ysize=0;
81
82 double sigma_x = 0.;
83 double sigma_y = 0.;
84
85 cpl_image *clean_flat;
86 cpl_image *eflat;
87 cpl_image *eflat_complex=NULL;
88 cpl_image *eflat_real=NULL;
89
90 cpl_image *filter_image=NULL;
91 cpl_image *filter_image_complex=NULL;
92 cpl_image *flat_real;
93
94
95 /*Clean flat from bad pixels if bpm exists*/
96 //clean_flat = cpl_image_duplicate(flat);
97
98 /*The algorithm uses float, so we have to cast here or change the algo...*/
99 cpl_type ima_type = cpl_image_get_type(ima);
100
101 clean_flat = cpl_image_cast(ima, CPL_TYPE_FLOAT);
102 cpl_detector_interpolate_rejected(clean_flat);
103
104
105 /* Expand the flat image using the mirror edges function */
106 eflat = hdrl_mirror_edges(clean_flat, mirrorx, mirrory);
107
108 if (clean_flat != NULL) {
109 cpl_image_delete(clean_flat);
110 clean_flat = NULL;
111 }
112
113 if(eflat == NULL){
114 cpl_msg_error(cpl_func,"Filter image is NULL");
115 return NULL;
116 }
117
118 xsize = cpl_image_get_size_x(eflat);
119 ysize = cpl_image_get_size_y(eflat);
120
121 sigma_x = gausfilt;
122 sigma_y = (double)(sigma_x * ysize) / xsize;
123
124
125 /* Generate a lowpass filter to be used in the FFT convolution */
126 filter_image = hdrl_gen_lowpass(xsize, ysize, sigma_x, sigma_y);
127 if(filter_image == NULL){
128 cpl_msg_error(cpl_func,"Filter image is NULL");
129 cpl_image_delete(eflat);
130 return NULL;
131 }
132
133 eflat_complex = cpl_image_new(xsize,ysize, CPL_TYPE_FLOAT_COMPLEX);
134 eflat_real = cpl_image_new(xsize,ysize, CPL_TYPE_FLOAT);
135 filter_image_complex =cpl_image_cast(filter_image,CPL_TYPE_FLOAT_COMPLEX);
136
137 /*Free memory*/
138 cpl_image_delete(filter_image);
139
140
141 /* Apply a forward FFT on the images */
142 cpl_fft_image(eflat_complex, eflat, CPL_FFT_FORWARD);
143 /*Free memory*/
144 cpl_image_delete(eflat);
145
146 /*Multiply the filter with the the FFT image */
147 cpl_image_multiply(eflat_complex,filter_image_complex);
148
149
150 /* Apply a backward FFT on the images */
151 cpl_fft_image(eflat_real, eflat_complex,CPL_FFT_BACKWARD);
152 /*Free memory*/
153 cpl_image_delete(eflat_complex);
154 cpl_image_delete(filter_image_complex);
155
156 /* Extract original image from the expanded image. */
157 flat_real = cpl_image_extract(eflat_real, mirrorx+1, mirrory+1,
158 xsize-mirrorx, ysize-mirrory);
159
160 if (flat_real == NULL) {
161 cpl_msg_error (cpl_func,"Real extracted image is NULL. <%s>", cpl_error_get_message());
162 return NULL;
163 }
164 cpl_image_delete(eflat_real);
165
166 cpl_image * out_double = cpl_image_cast(flat_real, ima_type);
167 cpl_image_delete(flat_real);
168 return out_double;
169}
170
171/*-------------------------------------------------------------------------*/
187/*--------------------------------------------------------------------------*/
188static cpl_image * hdrl_gen_lowpass(int xs, int ys, double sigma_x,
189 double sigma_y)
190{
191
192 int i= 0.0;
193 int j= 0.0;
194 int hlx= 0.0;
195 int hly = 0.0;
196 double x= 0.0;
197 double gaussval= 0.0;
198 float *data;
199
200 cpl_image *lowpass_image;
201
202
203 lowpass_image = cpl_image_new (xs, ys, CPL_TYPE_FLOAT);
204 if (lowpass_image == NULL) {
205 cpl_msg_error (cpl_func, "Cannot generate lowpass filter <%s>",cpl_error_get_message());
206 return NULL;
207 }
208
209 hlx = xs/2;
210 hly = ys/2;
211
212 data = cpl_image_get_data_float(lowpass_image);
213
214 /* Given an image with pixels 0<=i<N, 0<=j<M then the convolution image
215 has the following properties:
216
217 ima[0][0] = 1
218 ima[i][0] = ima[N-i][0] = exp (-0.5 * (i/sig_i)^2) 1<=i<N/2
219 ima[0][j] = ima[0][M-j] = exp (-0.5 * (j/sig_j)^2) 1<=j<M/2
220 ima[i][j] = ima[N-i][j] = ima[i][M-j] = ima[N-i][M-j]
221 = exp (-0.5 * ((i/sig_i)^2 + (j/sig_j)^2))
222 */
223
224 data[0] = (float)1.0;
225
226 /* first row */
227 for (i=1 ; i<=hlx ; i++) {
228 x = (double)i / sigma_x;
229 gaussval = (double)exp(-0.5*x*x);
230 data[i] = gaussval;
231 data[xs-i] = gaussval;
232 }
233
234 for (j=1; j<=hly ; j++) {
235 double y = (double)j / sigma_y;
236 /* first column */
237 data[j*xs] = (double)exp(-0.5*y*y);
238 data[(ys-j)*xs] = (double)exp(-0.5*y*y);
239
240 for (i=1 ; i<=hlx ; i++) {
241 /* Use internal symetries */
242 x = (double) i / sigma_x;
243 gaussval = (double)exp (-0.5*(x*x+y*y));
244 data[j*xs+i] = gaussval;
245 data[(j+1)*xs-i] = gaussval;
246 data[(ys-j)*xs+i] = gaussval;
247 data[(ys+1-j)*xs-i] = gaussval;
248
249 }
250 }
251
252 /* FIXME: for the moment, reset errno which is coming from exp()
253 in first for-loop at i=348. This is causing cfitsio to
254 fail when loading an extension image (bug in cfitsio too).
255 */
256 if(errno != 0)
257 errno = 0;
258
259 return lowpass_image;
260}
261
262
263
271/*--------------------------------------------------------------------------*/
272cpl_image *hdrl_mirror_edges(cpl_image *image, int dx, int dy)
273{
274
275 intptr_t xs = cpl_image_get_size_x(image);
276 intptr_t ys = cpl_image_get_size_y(image);
277 intptr_t xx = xs+(2*dx);
278 intptr_t yy = ys+(2*dy);
279 float * data = cpl_image_get_data_float(image);
280 cpl_image * big_image = cpl_image_new(xx, yy, CPL_TYPE_FLOAT);
281 float * out_data = cpl_image_get_data_float(big_image);
282
283 for (intptr_t j=0; j<ys ; j++){
284 intptr_t inrow = j*xs;
285 intptr_t outrow = (j+dy)*xx;
286
287 for (intptr_t i=0; i<xs ; i++){
288 out_data[outrow+dx+i] = data[inrow+i];
289 }
290
291 for (intptr_t i=0; i<dx; i++){
292 out_data[outrow+i] = data[inrow+dx-i-1];
293 out_data[outrow+xs+dx+i] = data[inrow+xs-i-1];
294
295 }
296 }
297
298 for (intptr_t j=0; j<dy ; j++) {
299
300 for (intptr_t i=0; i<xx; i++) {
301 out_data[j*xx+i] = out_data[(2*dy-j-1)*xx+i];
302 out_data[(yy-j-1)*xx+i] = out_data[(yy-2*dy+j)*xx+i];
303 }
304 }
305
306 return big_image;
307}
308
309
310
311cpl_image * hdrl_mime_image_polynomial_bkg(cpl_image * image,
312 int dim_X, int dim_Y,
313 cpl_matrix ** coeffs){
314 cpl_imagelist * imlist;
315 cpl_imagelist * bkg_imlist;
316 cpl_image * bkg_image;
317 cpl_image * bkg_image_origtype;
318
319 cpl_type ima_type;
320 if (image == NULL){
321 cpl_error_set_message(cpl_func,
322 CPL_ERROR_NULL_INPUT,
323 "Null input image provided");
324 return NULL;
325 }
326
327 ima_type = cpl_image_get_type(image);
328
329 imlist = cpl_imagelist_new();
330 bkg_imlist = cpl_imagelist_new();
331 cpl_imagelist_set(imlist, image, 0);
332
333 hdrl_mime_compute_polynomial_bkg(imlist, bkg_imlist, dim_X, dim_Y, coeffs);
334
335 cpl_imagelist_unwrap(imlist);
336 bkg_image = cpl_imagelist_unset(bkg_imlist, 0);
337 cpl_imagelist_delete(bkg_imlist);
338
339 bkg_image_origtype=cpl_image_cast(bkg_image, ima_type);
340 cpl_image_delete(bkg_image);
341
342 return bkg_image_origtype;
343
344}
345
346/*----------------------------------------------------------------------------*/
361/*----------------------------------------------------------------------------*/
363 const cpl_imagelist * images,
364 cpl_imagelist * bkg_images,
365 int dim_X,
366 int dim_Y,
367 cpl_matrix ** coeffs
368 )
369{
370
371 int n_images;
372 int n_x;
373 int n_y;
374 cpl_matrix * poly_tensors;
375 int n_tensor;
376 int im;
377 cpl_matrix * weights;
378
379 cpl_msg_debug(cpl_func, "Polynomial with X, Y dimensions %2d, %2d.",
380 dim_X, dim_Y);
381
382 /* Sanity Check of input data, and parameters */
383 cpl_error_ensure(images != NULL, CPL_ERROR_DATA_NOT_FOUND,
384 return CPL_ERROR_DATA_NOT_FOUND, "list of dithered images is empty");
385
386 cpl_error_ensure(cpl_imagelist_is_uniform(images) == 0, CPL_ERROR_INCOMPATIBLE_INPUT,
387 return CPL_ERROR_INCOMPATIBLE_INPUT, "input image list have non uniform data");
388
389 /*
390 Compute Dimensions
391 */
392 n_images = cpl_imagelist_get_size(images);
393 n_x = cpl_image_get_size_x(cpl_imagelist_get_const(images, 0));
394 n_y = cpl_image_get_size_y(cpl_imagelist_get_const(images, 0));
395
396 /*
397 Create tensor products of polynomials
398 */
399 poly_tensors = hdrl_mime_legendre_tensors_create(n_x, n_y, dim_X, dim_Y);
400 n_tensor = cpl_matrix_get_ncol(poly_tensors);
401 *coeffs = cpl_matrix_new(n_tensor, n_images);
402
403 weights = hdrl_mime_tensor_weights_create(n_x, n_y);
404
405 /*
406 Loop over each image to find the corresponding
407 sky background
408 */
409
410 for (im = 0; im < n_images; im++){
411 cpl_matrix * image_data;
412 cpl_matrix * bkg_image_data;
413 cpl_matrix * masked_tensors;
414 cpl_matrix * masked_image;
415 cpl_matrix * image_double_wrap;
416 cpl_matrix * coeff;
417 cpl_image * image;
418 cpl_image * bkg_image;
419 cpl_image * image_double;
420 cpl_image * bkg_image_float;
421 cpl_mask * mask_bin;
422 double alpha = 1.0e-10;
423
424 image_data = cpl_matrix_new(n_x*n_y, 1);
425 bkg_image_data = cpl_matrix_new(n_x*n_y, 1);
426 masked_image = cpl_matrix_new(n_x * n_y, 1);
427 masked_tensors = cpl_matrix_new(n_x * n_y, n_tensor);
428
429
430 /*
431 Load image, and mask
432 */
433
434 image = cpl_image_duplicate(cpl_imagelist_get_const(images, im));
435 mask_bin = cpl_image_get_bpm(image);
436 if (mask_bin == NULL) {
437 cpl_msg_info(cpl_func, "mask not available");
438 cpl_matrix_delete(poly_tensors);
439 cpl_matrix_delete(image_data);
440 cpl_matrix_delete(bkg_image_data);
441 cpl_matrix_delete(masked_image);
442 cpl_matrix_delete(masked_tensors);
443 cpl_image_delete(image);
444 return cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
445 }
446
447 image_double = cpl_image_cast(image, CPL_TYPE_DOUBLE);
448 image_double_wrap = cpl_matrix_wrap(n_x*n_y, 1,
449 cpl_image_get_data_double(image_double));
450 cpl_matrix_copy(image_data, image_double_wrap, 0, 0);
451
452 /*
453 Loop over all pixels, to reject from mask.
454 */
455 cpl_matrix_copy(masked_tensors, poly_tensors, 0, 0);
456 hdrl_mime_matrix_mask_rows(masked_tensors, mask_bin);
457 hdrl_mime_matrix_rescale_rows(masked_tensors, weights, masked_tensors);
458
459 cpl_matrix_copy(masked_image, image_data, 0, 0);
460 hdrl_mime_matrix_mask_rows(masked_image, mask_bin);
461 hdrl_mime_matrix_rescale_rows(masked_image, weights, masked_image);
462
463 /*
464 Find coefficients, agument the matrix of coefficients
465 */
466 coeff = hdrl_mime_linalg_solve_tikhonov(masked_tensors, masked_image,
467 alpha);
468 cpl_matrix_copy(*coeffs, coeff, 0, im);
469
470
471 /* hdrl_mime_matrix_formatted_dump(coeff, NULL); */
472
473 /*
474 Synthesize background, copy to image, and then to return imagelist
475 */
476 hdrl_mime_matrix_product(poly_tensors, coeff, bkg_image_data);
477 bkg_image = cpl_image_wrap_double(n_x, n_y,
478 cpl_matrix_get_data(bkg_image_data));
479 bkg_image_float = cpl_image_cast(bkg_image, CPL_TYPE_FLOAT);
480 cpl_imagelist_set(bkg_images, bkg_image_float, im);
481
482 cpl_matrix_delete(image_data);
483 cpl_matrix_delete(bkg_image_data);
484 cpl_matrix_delete(masked_image);
485 cpl_matrix_delete(masked_tensors);
486 cpl_matrix_delete(coeff);
487 cpl_image_delete(image);
488 cpl_image_delete(image_double);
489 cpl_matrix_unwrap(image_double_wrap);
490 cpl_image_unwrap(bkg_image);
491
492 }
493
494 /* save results, destroy whatever is not necessary, return */
495
496 cpl_matrix_delete(weights);
497 cpl_matrix_delete(poly_tensors);
498 return CPL_ERROR_NONE;
499}
500/*---------------------------------------------------------------------------*/
513/*---------------------------------------------------------------------------*/
514cpl_matrix *hdrl_mime_legendre_tensors_create(int nx, int ny, int npx, int npy)
515{
516 cpl_matrix *x;
517 cpl_matrix *y;
518 cpl_matrix *xpolys;
519 cpl_matrix *ypolys;
520 cpl_matrix *tensors;
521
522 double ax, ay, bx, by;
523
524/* testing input */
525 if (nx < 2 || ny < 2 || npx < 1 || npy < 1)
526 {
527 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
528 return NULL;
529 }
530
531/* setting parameters:
532 nx, ny numbers of the nodes in the x- and y-direction
533 ax, ay, bx, by endpoints of the intervals with nodes
534 */
535 ax = 0.0;
536 bx = nx - 1.0;
537
538 ay = 0.0;
539 by = ny - 1.0;
540
541/* creating equally spaced nodes */
542 x = hdrl_mime_matrix_linspace_create(nx, ax, bx);
543 y = hdrl_mime_matrix_linspace_create(ny, ay, by);
544
545/* creating the tensor products */
546 xpolys = hdrl_mime_legendre_polynomials_create(npx, ax, bx, x);
547 ypolys = hdrl_mime_legendre_polynomials_create(npy, ay, by, y);
549 xpolys);
550
551/* cleaning up */
552 cpl_matrix_delete(x);
553 cpl_matrix_delete(y);
554 cpl_matrix_delete(xpolys);
555 cpl_matrix_delete(ypolys);
556
557 return tensors;
558}
559/*---------------------------------------------------------------------------*/
571/*---------------------------------------------------------------------------*/
572cpl_matrix *hdrl_mime_matrix_linspace_create(int n, double a, double b)
573{
574 cpl_matrix *nodes;
575 double *data;
576 double h;
577 int i;
578
579/* testing input */
580 if (n < 2)
581 {
582 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
583 return NULL;
584 }
585
586/* allocating memory */
587 nodes = cpl_matrix_new(n, 1);
588 data = cpl_matrix_get_data(nodes);
589
590/* creating the nodes */
591 h = (b - a) / (n - 1);
592 for (i = 0; i < n; i++)
593 data[i] = a + i * h;
594 data[n - 1] = b;
595
596 return nodes;
597}
598/*---------------------------------------------------------------------------*/
616/*---------------------------------------------------------------------------*/
617cpl_matrix *hdrl_mime_legendre_polynomials_create(int npoly, double a, double b,
618 const cpl_matrix * x)
619{
620 cpl_matrix *polys;
621 double *m;
622 const double *mx;
623 double midpoint, scale;
624 int i, j, nr;
625
626/* testing input */
627 if (x == NULL)
628 {
629 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
630 return NULL;
631 }
632
633 if (npoly < 1 || a == b)
634 {
635 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
636 return NULL;
637 }
638
639/* The specific dimensions of the matrix x are not used, only its size. */
640 nr = cpl_matrix_get_nrow(x) * cpl_matrix_get_ncol(x);
641
642/* allocating memory */
643 polys = cpl_matrix_new(nr, npoly);
644
645/* initializing */
646 midpoint = 0.5 * (a + b);
647 scale = 2.0 / (b - a);
648
649/* filling the first column */
650 m = cpl_matrix_get_data(polys);
651 for (i = 0; i < nr; i++, m += npoly)
652 m[0] = 1.0;
653
654/* filling the second column */
655 m = cpl_matrix_get_data(polys);
656 mx = cpl_matrix_get_data_const(x);
657 if (npoly >= 2)
658 for (i = 0; i < nr; i++, m += npoly)
659 m[1] = scale * (mx[i] - midpoint);
660
661/* filling the remaining columns by recursion
662 * j P'_j(x) = (2j-1)xP'_{j-1}(x) - (j-1)P'_{j-2}(x)
663 * with P'(x) = P((2x - b - a)/(b - a)) for orthogonality in [-1, 1]->[a, b] */
664 m = cpl_matrix_get_data(polys);
665 for (i = 0; i < nr; i++, m += npoly)
666 {
667 double xi = scale * (mx[i] - midpoint);
668 for (j = 2; j < npoly; j++)
669 {
670 double alpha = (2.0 * j - 1.0) / j;
671 double beta = (j - 1.0) / j;
672 m[j] = alpha * xi * m[j - 1] - beta * m[j - 2];
673 }
674 }
675
676 return polys;
677}
678/*---------------------------------------------------------------------------*/
692/*---------------------------------------------------------------------------*/
694 cpl_matrix * mat1, const cpl_matrix * mat2)
695{
696 cpl_matrix *tensor;
697 cpl_matrix *repl1;
698 cpl_matrix *repl2;
699
700 int nc1, nc2;
701 int j1, j2, nc, col_count;
702
703/* testing input */
704 if (mat1 == NULL || mat2 == NULL)
705 {
706 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
707 return NULL;
708 }
709
710/* initializing */
711 nc1 = cpl_matrix_get_ncol(mat1);
712 nc2 = cpl_matrix_get_ncol(mat2);
713
714/* counting the admissible pairs */
715 nc = 0;
716 for (j1 = 0; j1 < nc1; j1++)
717 {
718 for (j2 = 0; j2 < nc2; j2++)
719 nc += (j1 * (nc2 - 1) + j2 * (nc1 - 1) <=
720 (nc1 - 1) * (nc2 - 1) ? 1 : 0);
721 }
722
723/* replicating the columns */
724 repl1 = cpl_matrix_new(cpl_matrix_get_nrow(mat1), nc);
725 repl2 = cpl_matrix_new(cpl_matrix_get_nrow(mat2), nc);
726 col_count = 0;
727
728 for (j1 = 0; j1 < nc1; j1++)
729 {
730 for (j2 = 0; j2 < nc2; j2++)
731 {
732 if (j1 * (nc2 - 1) + j2 * (nc1 - 1) <= (nc1 - 1) * (nc2 - 1))
733 {
734 hdrl_mime_matrix_copy_column(mat1, j1, repl1, col_count);
735 hdrl_mime_matrix_copy_column(mat2, j2, repl2, col_count);
736 col_count++;
737 }
738 }
739 }
740
741/* filling the matrix with the tensor products */
743 cpl_matrix_delete(repl1);
744 cpl_matrix_delete(repl2);
745
746 return tensor;
747}
748/*---------------------------------------------------------------------------*/
761/*---------------------------------------------------------------------------*/
762cpl_error_code hdrl_mime_matrix_copy_column(const cpl_matrix * mat1, int j_1,
763 cpl_matrix * mat2, int j_2)
764{
765 const double *m1;
766 double *m2;
767 int i, nr, nc1, nc2;
768
769/* testing input */
770 if (mat1 == NULL || mat2 == NULL)
771 return cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
772
773 if (cpl_matrix_get_nrow(mat1) != cpl_matrix_get_nrow(mat2))
774 return cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
775
776 if (j_1 < 0 || j_1 >= cpl_matrix_get_ncol(mat1) || j_2 < 0 ||
777 j_2 >= cpl_matrix_get_ncol(mat2))
778 return cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
779
780/* initializing */
781 nr = cpl_matrix_get_nrow(mat1);
782 nc1 = cpl_matrix_get_ncol(mat1);
783 nc2 = cpl_matrix_get_ncol(mat2);
784 m1 = cpl_matrix_get_data_const(mat1) + j_1;
785 m2 = cpl_matrix_get_data(mat2) + j_2;
786
787/* filling the column */
788 for (i = 0; i < nr; i++, m1 += nc1, m2 += nc2)
789 *m2 = *m1;
790
791 return CPL_ERROR_NONE;
792}
793/*---------------------------------------------------------------------------*/
805/*---------------------------------------------------------------------------*/
807 mat1, const cpl_matrix * mat2)
808{
809 cpl_matrix *tensor;
810
811 const double *m1;
812 double *t;
813 int nr1, nr2, nc;
814 int i1;
815
816/* testing input */
817 if (mat1 == NULL || mat2 == NULL)
818 {
819 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
820 return NULL;
821 }
822
823 if (cpl_matrix_get_ncol(mat1) != cpl_matrix_get_ncol(mat2))
824 {
825 cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
826 return NULL;
827 }
828
829/* initializing */
830 nr1 = cpl_matrix_get_nrow(mat1);
831 nr2 = cpl_matrix_get_nrow(mat2);
832 nc = cpl_matrix_get_ncol(mat1); /* nc = cpl_matrix_get_ncol(mat2); */
833
834/* allocating memory */
835 tensor = cpl_matrix_new(nr1 * nr2, nc);
836
837/* filling the matrix with the tensor products */
838 m1 = cpl_matrix_get_data_const(mat1);
839 t = cpl_matrix_get_data(tensor);
840 for (i1 = 0; i1 < nr1; i1++, m1 += nc)
841 {
842 int i2;
843 const double *m2;
844 m2 = cpl_matrix_get_data_const(mat2);
845 for (i2 = 0; i2 < nr2; i2++, m2 += nc, t += nc)
846 {
847 int j;
848 for (j = 0; j < nc; j++)
849 t[j] = m1[j] * m2[j];
850 }
851 }
852
853 return tensor;
854}
855
856/*---------------------------------------------------------------------------*/
867/*---------------------------------------------------------------------------*/
868cpl_matrix *hdrl_mime_tensor_weights_create(int nx, int ny)
869{
870 cpl_matrix *x;
871 cpl_matrix *y;
872 cpl_matrix *weights;
873
874 double *m;
875 double ax, ay, bx, by, v;
876 int i;
877
878/* testing input */
879 if (nx < 2 || ny < 2)
880 {
881 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
882 return NULL;
883 }
884
885/* setting parameters:
886 nx, ny numbers of the nodes in the x- and y-direction
887 */
888 bx = 1.0 - 1.0 / nx;
889 ax = -bx;
890
891 by = 1.0 - 1.0 / ny;
892 ay = -by;
893
894/* creating equally spaced nodes */
895 x = hdrl_mime_matrix_linspace_create(nx, ax, bx);
896 y = hdrl_mime_matrix_linspace_create(ny, ay, by);
897
898/* creating the weights as the tensor product */
899 m = cpl_matrix_get_data(x);
900 for (i = 0; i < nx; i++)
901 {
902 v = m[i];
903 v = 1.0 / sqrt(1.0 - v * v);
904 m[i] = sqrt(v);
905 }
906
907 m = cpl_matrix_get_data(y);
908 for (i = 0; i < ny; i++)
909 {
910 v = m[i];
911 v = 1.0 / sqrt(1.0 - v * v);
912 m[i] = sqrt(v);
913 }
914
915/* switch
916 no weights: 1
917 with weights: 0 */
918 if (1)
919 {
920 cpl_matrix_fill(x, 1.0);
921 cpl_matrix_fill(y, 1.0);
922 }
923
925
926/* cleaning up */
927 cpl_matrix_delete(x);
928 cpl_matrix_delete(y);
929
930 return weights;
931}
932/*---------------------------------------------------------------------------*/
944/*---------------------------------------------------------------------------*/
945cpl_error_code hdrl_mime_matrix_mask_rows(cpl_matrix * mat, const cpl_mask * mask)
946{
947 double *m;
948 const cpl_binary *mk;
949 int i, j, nr, nc;
950
951/* testing input */
952 if (mat == NULL || mask == NULL)
953 return cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
954
955 if (cpl_matrix_get_nrow(mat) !=
956 cpl_mask_get_size_x(mask) * cpl_mask_get_size_y(mask))
957 return cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
958
959/* initializing */
960 nr = cpl_matrix_get_nrow(mat);
961 nc = cpl_matrix_get_ncol(mat);
962 m = cpl_matrix_get_data(mat);
963 mk = cpl_mask_get_data_const(mask);
964
965/* updating the rows */
966 for (i = 0; i < nr; i++, mk++, m += nc)
967 {
968 if (*mk == CPL_BINARY_1)
969 {
970 for (j = 0; j < nc; j++)
971 m[j] = 0.0;
972 }
973 }
974
975 return CPL_ERROR_NONE;
976}
977/*---------------------------------------------------------------------------*/
990/*---------------------------------------------------------------------------*/
991cpl_error_code hdrl_mime_matrix_rescale_rows(const cpl_matrix * mat,
992 const cpl_matrix * d, cpl_matrix * dmat)
993{
994 const double *m;
995 const double *di;
996 double *dm;
997 int i, j, nr, nc;
998
999/* testing input */
1000 if (mat == NULL || d == NULL || dmat == NULL)
1001 {
1002 return cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1003 }
1004
1005 if (cpl_matrix_get_nrow(mat) !=
1006 cpl_matrix_get_nrow(d) * cpl_matrix_get_ncol(d))
1007 {
1008 return cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1009 }
1010
1011 if (cpl_matrix_get_ncol(mat) != cpl_matrix_get_ncol(dmat) ||
1012 cpl_matrix_get_nrow(mat) != cpl_matrix_get_nrow(dmat))
1013 {
1014 return cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1015 }
1016
1017/* initializing */
1018 nr = cpl_matrix_get_nrow(mat);
1019 nc = cpl_matrix_get_ncol(mat);
1020
1021 m = cpl_matrix_get_data_const(mat);
1022 di = cpl_matrix_get_data_const(d);
1023 dm = cpl_matrix_get_data(dmat);
1024
1025/* multiplying the rows by the d[i]'s */
1026 for (i = 0; i < nr; i++, m += nc, dm += nc)
1027 {
1028 for (j = 0; j < nc; j++)
1029 dm[j] = di[i] * m[j];
1030 }
1031
1032 return CPL_ERROR_NONE;
1033}
1034/*---------------------------------------------------------------------------*/
1050/*---------------------------------------------------------------------------*/
1051cpl_matrix *hdrl_mime_linalg_solve_tikhonov(const cpl_matrix * mat,
1052 const cpl_matrix * rhs, double alpha)
1053{
1054 cpl_matrix *normal;
1055 cpl_matrix *solution;
1056 cpl_error_code error;
1057
1058/* testing input */
1059 if (mat == NULL || rhs == NULL)
1060 {
1061 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1062 return NULL;
1063 }
1064
1065 if (cpl_matrix_get_nrow(mat) != cpl_matrix_get_nrow(rhs))
1066 {
1067 cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1068 return NULL;
1069 }
1070
1071/* creating the normal equations and computing the Cholesky decomposition */
1072 normal = hdrl_mime_linalg_normal_equations_create(mat, alpha);
1073 error = cpl_matrix_decomp_chol(normal);
1074
1075 if (error != CPL_ERROR_NONE)
1076 {
1077 cpl_matrix_delete(normal);
1078 return NULL;
1079 }
1080
1081/* solving the normal equations and cleaning up */
1083 error = cpl_matrix_solve_chol(normal, solution);
1084 cpl_matrix_delete(normal);
1085
1086 if (error != CPL_ERROR_NONE)
1087 {
1088 cpl_matrix_delete(solution);
1089 return NULL;
1090 }
1091
1092 return solution;
1093}
1094/*---------------------------------------------------------------------------*/
1106/*---------------------------------------------------------------------------*/
1107cpl_matrix *hdrl_mime_linalg_normal_equations_create(const cpl_matrix * mat,
1108 double alpha)
1109{
1110 cpl_matrix *normal;
1111 const double *m;
1112 double *p;
1113 double sum;
1114 int nr, nc;
1115 int i, j, k;
1116
1117/* testing input */
1118 if (mat == NULL)
1119 {
1120 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1121 return NULL;
1122 }
1123
1124 if (alpha < 0.0)
1125 {
1126 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
1127 return NULL;
1128 }
1129
1130/* initializing */
1131 nr = cpl_matrix_get_nrow(mat);
1132 nc = cpl_matrix_get_ncol(mat);
1133
1134/* allocating memory */
1135 normal = cpl_matrix_new(nc, nc);
1136
1137/* filling the normal matrix */
1138 p = cpl_matrix_get_data(normal);
1139 for (i = 0; i < nc; i++, p += nc)
1140 {
1141 for (j = i; j < nc; j++)
1142 {
1143 m = cpl_matrix_get_data_const(mat);
1144 sum = 0.0;
1145 for (k = 0; k < nr; k++, m += nc)
1146 sum += m[i] * m[j];
1147 p[j] = sum;
1148 }
1149 }
1150
1151/* updating the diagonal */
1152 p = cpl_matrix_get_data(normal);
1153 for (i = 0; i < nc; i++)
1154 p[nc * i + i] += alpha;
1155
1156 return normal;
1157}
1158/*---------------------------------------------------------------------------*/
1172/*---------------------------------------------------------------------------*/
1174 mat1, const cpl_matrix * mat2)
1175{
1176 cpl_matrix *product;
1177
1178 const double *m1;
1179 const double *m2;
1180 double *p;
1181 double sum;
1182
1183 int nr, nc, common;
1184 int i, j, k;
1185
1186/* testing input */
1187 if (mat1 == NULL || mat2 == NULL)
1188 {
1189 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1190 return NULL;
1191 }
1192
1193 if (cpl_matrix_get_nrow(mat1) != cpl_matrix_get_nrow(mat2))
1194 {
1195 cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1196 return NULL;
1197 }
1198
1199/* initializing */
1200 nr = cpl_matrix_get_ncol(mat1);
1201 nc = cpl_matrix_get_ncol(mat2);
1202 common = cpl_matrix_get_nrow(mat1);
1203 /* common = cpl_matrix_get_nrow(mat2); */
1204
1205/* allocating memory */
1206 product = cpl_matrix_new(nr, nc);
1207 p = cpl_matrix_get_data(product);
1208
1209/* filling the product matrix */
1210 for (i = 0; i < nr; i++, p += nc)
1211 {
1212 for (j = 0; j < nc; j++)
1213 {
1214 m1 = cpl_matrix_get_data_const(mat1);
1215 m2 = cpl_matrix_get_data_const(mat2);
1216 sum = 0.0;
1217 for (k = 0; k < common; k++, m1 += nr, m2 += nc)
1218 sum += m1[i] * m2[j];
1219 p[j] = sum;
1220 }
1221 }
1222
1223 return product;
1224}
1225/*---------------------------------------------------------------------------*/
1240/*---------------------------------------------------------------------------*/
1241cpl_error_code hdrl_mime_matrix_product(const cpl_matrix * mat1,
1242 const cpl_matrix * mat2, cpl_matrix * product)
1243{
1244 const double *m1;
1245 const double *m2;
1246 double *p;
1247 double sum;
1248
1249 int nr, nc, common;
1250 int i, j, k;
1251
1252/* testing input */
1253 if (mat1 == NULL || mat2 == NULL || product == NULL)
1254 {
1255 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1256 return CPL_ERROR_NONE;
1257 }
1258
1259 if (cpl_matrix_get_ncol(mat1) != cpl_matrix_get_nrow(mat2) ||
1260 cpl_matrix_get_nrow(mat1) != cpl_matrix_get_nrow(product) ||
1261 cpl_matrix_get_ncol(mat2) != cpl_matrix_get_ncol(product))
1262 {
1263 cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1264 return CPL_ERROR_NONE;
1265 }
1266
1267/* initializing */
1268 nr = cpl_matrix_get_nrow(mat1);
1269 nc = cpl_matrix_get_ncol(mat2);
1270 common = cpl_matrix_get_ncol(mat1); /* common = cpl_matrix_get_nrow(mat2); */
1271
1272/* filling the product matrix */
1273 m1 = cpl_matrix_get_data_const(mat1);
1274 p = cpl_matrix_get_data(product);
1275 for (i = 0; i < nr; i++, m1 += cpl_matrix_get_ncol(mat1), p += nc)
1276 {
1277 for (j = 0; j < nc; j++)
1278 {
1279 m2 = cpl_matrix_get_data_const(mat2);
1280 sum = 0.0;
1281 for (k = 0; k < common; k++, m2 += cpl_matrix_get_ncol(mat2))
1282 sum += m1[k] * m2[j];
1283 p[j] = sum;
1284 }
1285 }
1286
1287 return CPL_ERROR_NONE;
1288}
1289
1290
1291
cpl_error_code hdrl_mime_compute_polynomial_bkg(const cpl_imagelist *images, cpl_imagelist *bkg_images, int dim_X, int dim_Y, cpl_matrix **coeffs)
Fit smooth background for a list of images.
cpl_matrix * hdrl_mime_legendre_tensors_create(int nx, int ny, int npx, int npy)
Create tensor products of Legendre polynomials.
cpl_error_code hdrl_mime_matrix_mask_rows(cpl_matrix *mat, const cpl_mask *mask)
Fill matrix rows with zeros as indicated by a mask.
cpl_matrix * hdrl_mime_linalg_solve_tikhonov(const cpl_matrix *mat, const cpl_matrix *rhs, double alpha)
Solve an overdetermined linear system in the least-squares sense.
cpl_error_code hdrl_mime_matrix_copy_column(const cpl_matrix *mat1, int j_1, cpl_matrix *mat2, int j_2)
Copy a column from one matrix to another matrix.
cpl_matrix * hdrl_mime_linalg_pairwise_column_tensor_products_create(const cpl_matrix *mat1, const cpl_matrix *mat2)
Create selected pairwise tensor products of the columns of two matrices.
cpl_matrix * hdrl_mime_matrix_product_left_transpose_create(const cpl_matrix *mat1, const cpl_matrix *mat2)
Create the product of the transpose of a matrix with another matrix.
cpl_error_code hdrl_mime_matrix_rescale_rows(const cpl_matrix *mat, const cpl_matrix *d, cpl_matrix *dmat)
Multiply the rows of a matrix by given factors.
cpl_matrix * hdrl_mime_tensor_weights_create(int nx, int ny)
Create tensor product weights.
cpl_matrix * hdrl_mime_legendre_polynomials_create(int npoly, double a, double b, const cpl_matrix *x)
Create the Legendre polynomial basis on the interval (a,b).
cpl_image * hdrl_get_spatial_freq(cpl_image *ima, double gausfilt, int mirrorx, int mirrory)
Get low spatial frequency componenets from an image using the FFTW.
cpl_matrix * hdrl_mime_matrix_linspace_create(int n, double a, double b)
Create equally spaced nodes.
cpl_matrix * hdrl_mime_linalg_normal_equations_create(const cpl_matrix *mat, double alpha)
Create the matrix transpose(A) * A + alpha for given A and alpha.
cpl_matrix * hdrl_mime_linalg_tensor_products_columns_create(const cpl_matrix *mat1, const cpl_matrix *mat2)
Create the tensor products of the columns of two matrices.
cpl_error_code hdrl_mime_matrix_product(const cpl_matrix *mat1, const cpl_matrix *mat2, cpl_matrix *product)
Fill a matrix with the product of two given matrices.