GRAVI Pipeline Reference Manual 1.9.6
Loading...
Searching...
No Matches
gravi_acqcam.c
Go to the documentation of this file.
1/* $Id: gravi_vis.c,v 1.10 2014/11/12 15:10:40 nazouaoui Exp $
2 *
3 * This file is part of the GRAVI Pipeline
4 * Copyright (C) 2002,2003 European Southern Observatory
5 *
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 2 of the License, or
9 * (at your option) any later version.
10 *
11 * This program is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 */
20
33/*
34 * History :
35 * 26/11/2018 Read parameter roof_x, roof_y, sot_x, spot_y and roof_pos from calibration file
36 * 28/11/2018 Read parameter plate scale from calibration file
37 * 03/12/2018 Read gravi_acqcam_pupil parameter from calibration file
38 * 10/01/2019 Fix compiler Warnings : inititialize with NULL : roof_x, roof_y, sot_x, spot_y and roof_pos
39 * unused parameter nsy
40 * 12/09/2019 EKW Flip parameter from cpl_vector_set (output, GRAVI_SPOT_DIODE+0, 0.122);
41 * cpl_vector_set (output, GRAVI_SPOT_DIODE+1, 0.158)
42 * to cpl_vector_set (output, GRAVI_SPOT_DIODE+0, 0.158);
43 * cpl_vector_set (output, GRAVI_SPOT_DIODE+1, 0.122)
44 * 18/09/2019 EkW Double spot fitting as Ollis analysis showed problems
45 * 16/10/2019 EkW Remove debug statements
46 * 06/11/2019 EKW Add Franks double pscale = sep ? rho_in/sep : scale (instead 0) ;
47 * 06/11/2019 EKW Commit with right PIPE-8675 as despite the line abovem, all changes are related to it
48 * 15/07/2020 EKW JIRA PIPE-9123 : Pipeline faile with Strehl=NaN
49 */
50/*-----------------------------------------------------------------------------
51 Includes
52 -----------------------------------------------------------------------------*/
53
54#ifdef HAVE_CONFIG_H
55#include <config.h>
56#endif
57
58#include <cpl.h>
59#include <string.h>
60#include <stdio.h>
61#include <math.h>
62#include <time.h>
63
64#include <hdrl_strehl.h>
65
66#include "gravi_dfs.h"
67#include "gravi_data.h"
68#include "gravi_pfits.h"
69
70#include "gravi_cpl.h"
71#include "gravi_utils.h"
72
73#include "gravi_acqcam.h"
74
75
76/*-----------------------------------------------------------------------------
77 Defines & Macros
78 -----------------------------------------------------------------------------*/
79
80/* Number of parameters in the model 'gravi_acqcam_spot'
81 * And position of parameters */
82#define GRAVI_SPOT_NTEL 4
83#define GRAVI_SPOT_NLENS 4
84#define GRAVI_SPOT_NSPOT 4
85#define GRAVI_SPOT_NFOCUS 21
86#define GRAVI_SPOT_NSEARCH 45
87#define GRAVI_SPOT_SWINDOW 28
88
89#define GRAVI_SPOT_NA 30
90#define GRAVI_SPOT_SUB 0
91#define GRAVI_SPOT_ANGLE 8
92#define GRAVI_SPOT_SCALE 9
93#define GRAVI_SPOT_DIODE 10
94#define GRAVI_SPOT_FWHM 13
95#define GRAVI_SPOT_FLUX 14
96
97#define GRAVI_ACQ_PUP_FLUX 1e6
98
99#define CHECK_BIT(var,pos) ((var) & (1<<(pos)))
100
101/*-----------------------------------------------------------------------------
102 Private prototypes
103 -----------------------------------------------------------------------------*/
104
105double exp1 (double x);
106double sin1 (double x);
107int gravi_acqcam_xy_diode (const double v[], double *xd, double *yd);
108
109//static int gravi_acqcam_spot (const double x_in[], const double v[], double *result);
110//static int gravi_acqcam_xy_sub (const double v[], double *xsub, double *ysub);
111
112double gravi_acqcam_z2meter (double PositionPixels, gravi_data *static_param_data);
113
114cpl_error_code gravi_acqcam_pupil (cpl_image * mean_img,
115 cpl_imagelist * acqcam_imglist,
116 cpl_propertylist * header,
117 cpl_table * acqcam_table,
118 cpl_propertylist * o_header,
119 gravi_data *static_param_data);
120
121cpl_error_code gravi_acqcam_field (cpl_image * mean_img,
122 cpl_imagelist * acqcam_imglist,
123 cpl_propertylist * header,
124 cpl_table * acqcam_table,
125 cpl_propertylist * o_header,
126 gravi_data *static_param_data);
127
128cpl_error_code gravi_acq_fit_gaussian (cpl_image * img, double *x, double *y,
129 double *ex, double *ey, cpl_size size);
130
131cpl_error_code gravi_acq_measure_strehl(cpl_image * img, double x, double y,
132 double pscale, double *SR, cpl_propertylist * header);
133
134cpl_error_code gravi_acq_measure_max(cpl_image * img, double x, double y, double size, double * img_max);
135
136cpl_error_code gravi_image_fft_correlate (cpl_image *ia, cpl_image *ib, cpl_size *xd, cpl_size *yd);
137
138cpl_error_code gravi_acqcam_pupil_v2 (cpl_image * mean_img,
139 cpl_imagelist * acqcam_imglist,
140 cpl_propertylist * header,
141 cpl_table * acqcam_table,
142 cpl_propertylist * o_header,
143 gravi_data *static_param_data);
144
145cpl_error_code gravi_acqcam_clean_pupil_v2(cpl_imagelist * acqcam_imglist, cpl_imagelist * pupilImage_filtered, const cpl_size ury);
146
147cpl_error_code gravi_acqcam_select_good_frames_v2(cpl_imagelist * acqcam_imglist, cpl_imagelist * pupilImage_onFrames, cpl_array * good_frames);
148
149cpl_error_code gravi_acqcam_get_pup_ref_v2 (cpl_propertylist * header,
150 cpl_bivector * diode_pos_subwindow);
151
152cpl_error_code gravi_acqcam_get_diode_ref_v2 (cpl_propertylist * header,
153 cpl_array * good_frames,
154 cpl_vector * scale_vector,
155 cpl_bivector ** diode_pos_telescope,
156 int nrow_on);
157
158cpl_error_code gravi_acqcam_get_diode_theoretical_v2(cpl_bivector * diode_pos_subwindow,
159 cpl_bivector ** diode_pos_telescope,
160 cpl_bivector ** diode_pos_theoretical,
161 cpl_size nrow_on, int ury);
162
163cpl_error_code gravi_acqcam_spot_imprint_v2(cpl_image *mean_img, cpl_bivector ** diode_pos_offset, cpl_bivector ** diode_pos_theoretical, int ury);
164
165
166cpl_error_code gravi_acqcam_perform_shiftandadd_v2(cpl_imagelist * pupilImage_onFrames, cpl_imagelist ** pupilImage_shiftandadd, cpl_array * good_frames,
167 cpl_vector * focus_value,
168 cpl_bivector ** diode_pos_theoretical,
169 cpl_bivector ** diode_pos_offset, cpl_size nrow_on);
170
171cpl_error_code gravi_acqcam_get_pupil_offset_v2(cpl_imagelist ** pupilImage_shiftandadd,
172 cpl_array * good_frames,
173 cpl_bivector ** diode_pos_offset, cpl_propertylist * o_header,
174 cpl_size nrow_on);
175
176cpl_error_code gravi_acqcam_set_pupil_table_v2(cpl_table * acqcam_table, cpl_propertylist * header, cpl_propertylist * o_header, cpl_array * good_frames,
177 cpl_vector* scale_vector, cpl_array * bad_frames_short, cpl_bivector ** diode_pos_offset , cpl_vector * focus_value, gravi_data *static_param_data);
178
179cpl_image * gravi_image_extract(cpl_image * image_in, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury);
180
181double gravi_acqcam_defocus_scaling(int focus);
182
183/* This global variable optimises the computation
184 * of partial derivative on fitted parameters */
185const int * GRAVI_LVMQ_FREE;
186const int * GRAVI_LVMQ_FREE = NULL;
187
188
189/*-----------------------------------------------------------------------------
190 Functions code
191 -----------------------------------------------------------------------------*/
192
193/*----------------------------------------------------------------------------*/
206/*----------------------------------------------------------------------------*/
207
208cpl_error_code gravi_preproc_acqcam (gravi_data *output_data,
209 gravi_data *input_data,
210 gravi_data *bad_map)
211{
213 cpl_ensure_code (output_data, CPL_ERROR_NULL_INPUT);
214 cpl_ensure_code (input_data, CPL_ERROR_NULL_INPUT);
215
216 /* Check if extension exist */
218 gravi_msg_warning (cpl_func,"Cannot preproc the ACQCAM, no data in file");
219 return CPL_ERROR_NONE;
220 }
222 gravi_msg_warning (cpl_func,"Cannot preproc the ACQCAM, no badpixel in BAD");
223 return CPL_ERROR_NONE;
224 }
225
226 /* get image of badpixels */
227 cpl_image * badpix_img = gravi_data_get_img (bad_map, GRAVI_IMAGING_DATA_ACQ_EXT);
228 CPLCHECK_MSG ("Cannot get BAD map for ACQ");
229
230 /* Get the imagelist */
231 cpl_imagelist * imglist;
232 imglist = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
233 CPLCHECK_MSG ("Cannot get image for ACQ");
234
235 /* check for image size */
236 cpl_image * image0 = cpl_imagelist_get (imglist, 0);
237 if ( (cpl_image_get_size_x (image0) != cpl_image_get_size_x (badpix_img) ) ||
238 (cpl_image_get_size_y (image0) != cpl_image_get_size_y (badpix_img) ) ) {
239 gravi_msg_warning (cpl_func,"Cannot preproc the ACQCAM. Bad pixel mask does not have correct size.");
240 return CPL_ERROR_NONE;
241 }
242
243 /* Construct a mask of badpixels */
244 cpl_mask * badpix_mask = cpl_mask_threshold_image_create (badpix_img, 0.5, 10000);
245
246 /* Allocate new memory */
247 imglist = cpl_imagelist_duplicate (imglist);
248
249
250 /*
251 * Loop on images to cleanup-badpixels
252 */
253 cpl_size nrow = cpl_imagelist_get_size (imglist);
254 for (cpl_size row = 0; row < nrow; row++) {
255
256 /* Get image */
257 cpl_image * img = cpl_imagelist_get (imglist, row);
258
259 /* Cleanup-badpixel */
260 cpl_image_reject_from_mask (img, badpix_mask);
261 cpl_detector_interpolate_rejected (img);
262 CPLCHECK_MSG ("Cannot clean-up badpixel");
263 }
264
265 FREE (cpl_mask_delete, badpix_mask);
266
267 gravi_data_add_cube (output_data, NULL, GRAVI_IMAGING_DATA_ACQ_EXT, imglist);
268
270 return CPL_ERROR_NONE;
271}
272
273/*----------------------------------------------------------------------------*/
274
275/* Fast sin function */
276double sin1 (double x)
277{
278 /* Within -pi +pi*/
279 while (x >= CPL_MATH_PI) x -= CPL_MATH_2PI;
280 while (x <= -CPL_MATH_PI) x += CPL_MATH_2PI;
281
282 /* Within -pi/2 +pi/2*/
283 if (x > CPL_MATH_PI_2) x = CPL_MATH_PI - x;
284 if (x < -CPL_MATH_PI_2) x = -CPL_MATH_PI - x;
285
286 double x2 = x*x, x3 = x2*x, x5 = x3*x2;
287 return 0.9996949 * x - 0.1656700 * x3 + 0.0075134 * x5;
288}
289
290/* Fast exp function */
291double exp1 (double x) {
292 x = 1.0 + x / 256.0;
293 x *= x; x *= x; x *= x; x *= x;
294 x *= x; x *= x; x *= x; x *= x;
295 return x;
296}
297
298/* Compute the 4 diode positions from {angle, scaling, dx, dy}
299 * The model assume the 4 diodes form a rectangle centered on
300 * the pupil. Hence this model has a 180deg degeneracy */
301int gravi_acqcam_xy_diode (const double v[], double *xd, double *yd)
302{
303 double dx = v[GRAVI_SPOT_DIODE+0];
304 double dy = v[GRAVI_SPOT_DIODE+1];
305
306 /* Angle */
307 double ang = (v[GRAVI_SPOT_ANGLE] - v[GRAVI_SPOT_DIODE+2])* CPL_MATH_RAD_DEG;
308 double sang = sin1 (ang) * v[GRAVI_SPOT_SCALE];
309 double cang = sin1 (ang + CPL_MATH_PI_2) * v[GRAVI_SPOT_SCALE];
310
311 /* Diode arrangement */
312 xd[0] = -cang * dx + sang * dy;
313 xd[1] = -cang * dx - sang * dy;
314 xd[2] = cang * dx - sang * dy;
315 xd[3] = cang * dx + sang * dy;
316
317 yd[0] = cang * dy + sang * dx;
318 yd[1] = -cang * dy + sang * dx;
319 yd[2] = -cang * dy - sang * dx;
320 yd[3] = cang * dy - sang * dx;
321
322 return 0;
323}
324
325#if 0
326/* Compute the 4 sub-aperture positions from the sub-aperture modes */
327static int gravi_acqcam_xy_sub (const double v[], double *xsub, double *ysub)
328{
329 /* Sub-apperture arrangement */
330 const double * vd = v + GRAVI_SPOT_SUB;
331 xsub[0] = vd[0] + vd[1] + vd[2] + vd[3];
332 xsub[1] = vd[0] - vd[1] + vd[2] - vd[3];
333 xsub[2] = vd[0] + vd[1] - vd[2] - vd[3];
334 xsub[3] = vd[0] - vd[1] - vd[2] + vd[3];
335
336 ysub[0] = vd[4] + vd[5] + vd[6] + vd[7];
337 ysub[1] = vd[4] - vd[5] + vd[6] - vd[7];
338 ysub[2] = vd[4] + vd[5] - vd[6] - vd[7];
339 ysub[3] = vd[4] - vd[5] - vd[6] + vd[7];
340
341 return 0;
342}
343
344/*----------------------------------------------------------------------------*/
345
346static int gravi_acqcam_spot (const double x_in[], const double v[], double *result)
347{
348 *result = 0.0;
349
350 /* Static parameters */
351 double xd[4], yd[4], xsub[4], ysub[4];
352 gravi_acqcam_xy_diode (v, xd, yd);
353 gravi_acqcam_xy_sub (v, xsub, ysub);
354
355 double fwhm2 = v[GRAVI_SPOT_FWHM];
356
357 /* Loop on diode and appertures.
358 * The capture range is 2.FWHM */
359 for (int diode = 0; diode < 4 ; diode++) {
360 for (int sub = 0; sub < 4 ; sub++) {
361 double xf = (x_in[0] - xsub[sub] - xd[diode]);
362 double yf = (x_in[1] - ysub[sub] - yd[diode]);
363 double dist = (xf*xf + yf*yf) / fwhm2;
364 if (dist < 4.0) *result += v[GRAVI_SPOT_FLUX+sub*4+diode] * exp1 (-dist);
365 }
366 }
367
368 return 0;
369}
370
371/*----------------------------------------------------------------------------*/
372static int gravi_acqcam_spot_dfda (const double x_in[], const double v[], double result[])
373{
374 double next = 0.0, here = 0.0, epsilon = 1e-8;
375
376 double vlocal[GRAVI_SPOT_NA];
377 memcpy (vlocal, v, sizeof(double)*GRAVI_SPOT_NA);
378
379 /* Compute value in-place */
380 gravi_acqcam_spot (x_in, vlocal, &here);
381
382 /* Fill with zeros */
383 for (int a = 0; a < GRAVI_SPOT_NA; a++) result[a] = 0.0;
384
385 /* Loop on parameters to compute finite differences
386 * FIXME: The analytical derivative may be faster, but
387 * wasn't true in first tests, thus keep these. */
388 for (int a = 0; a < GRAVI_SPOT_FLUX; a++) {
389 if (GRAVI_LVMQ_FREE[a] != 0) {
390 vlocal[a] += epsilon;
391 gravi_acqcam_spot (x_in, vlocal, &next);
392 vlocal[a] -= epsilon;
393 result[a] = (next - here) / epsilon;
394 }
395 }
396
397 /* The intensities are trivial analytic derivative */
398 for (int a = GRAVI_SPOT_FLUX; a < GRAVI_SPOT_NA; a++) {
399 if (GRAVI_LVMQ_FREE[a] != 0) {
400 result[a] = here / v[a];
401 }
402 }
403
404 return 0;
405}
406#endif
407
408/*----------------------------------------------------------------------------*/
421/*----------------------------------------------------------------------------*/
422
423cpl_error_code gravi_acq_measure_strehl(cpl_image * img, double x, double y,
424 double pscale, double *SR,
425 cpl_propertylist * header)
426{
428 cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
429 cpl_ensure_code (SR, CPL_ERROR_NULL_INPUT);
430 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
431 cpl_ensure_code (pscale > 0, CPL_ERROR_ILLEGAL_INPUT);
432
433 cpl_size hw = 45;
434 cpl_ensure_code (x > hw, CPL_ERROR_ILLEGAL_INPUT);
435 cpl_ensure_code (y > hw, CPL_ERROR_ILLEGAL_INPUT);
436
437 /* Structure to define the parameters
438 for fitting strehl in image */
439 hdrl_parameter * strehl_params;
440
441 /* Hardcoded Mirror diameters */
442 const char * telname = gravi_conf_get_telname (0, header);
443 CPLCHECK ("Cannot get telescope name");
444
445 if (telname[0] == 'U') {
446 strehl_params =
447 hdrl_strehl_parameter_create (1.8e-6, 8.0/2.0, 1.126/2.0,
448 pscale*1e-3, pscale*1e-3,
449 0.8, 0.8, 1.0);
450 } else if (telname[0] == 'A') {
451 strehl_params =
452 hdrl_strehl_parameter_create (1.8e-6, 1.8/2.0, 0.138/2.0,
453 pscale*1e-3, pscale*1e-3,
454 0.8, 0.8, 1.0);
455 } else {
456 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
457 "Cannot recognise the telescope");
458 return CPL_ERROR_ILLEGAL_INPUT;
459 }
460
461 /* Extract sub-image */
462 cpl_image * sub_image = cpl_image_extract (img, x-hw, y-hw, x+hw, y+hw);
463 hdrl_image * sub_hdrl = hdrl_image_create (sub_image, NULL);
464
465 /* Run Strehl algorithm from HDRL package */
466 hdrl_strehl_result strehl;
467 strehl = hdrl_strehl_compute (sub_hdrl, strehl_params);
468 *SR = (double) strehl.strehl_value.data;
469
470 /* Delete */
471 FREE (cpl_image_delete, sub_image);
472 FREE (hdrl_image_delete, sub_hdrl);
473 FREE (hdrl_parameter_delete, strehl_params);
474
476 return CPL_ERROR_NONE;
477}
478
479/*----------------------------------------------------------------------------*/
480/*
481 * @brief measure maximum of the source at the given location
482 */
483/*----------------------------------------------------------------------------*/
484
485cpl_error_code gravi_acq_measure_max(cpl_image * img, double x, double y, double size, double * img_max)
486{
488 cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
489 cpl_ensure_code (img_max, CPL_ERROR_NULL_INPUT);
490
491 cpl_image * small_img = cpl_image_extract(img, x-size, y-size, x+size, y+size);
492 *img_max = cpl_image_get_max(small_img);
493 cpl_image_delete(small_img);
494
496 return CPL_ERROR_NONE;
497}
498
499
500/*----------------------------------------------------------------------------*/
514/*----------------------------------------------------------------------------*/
515
516cpl_error_code gravi_acq_fit_gaussian (cpl_image * img, double *x, double *y,
517 double *ex, double *ey, cpl_size size)
518{
520 cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
521 cpl_ensure_code (x, CPL_ERROR_NULL_INPUT);
522 cpl_ensure_code (y, CPL_ERROR_NULL_INPUT);
523 cpl_ensure_code (ex, CPL_ERROR_NULL_INPUT);
524 cpl_ensure_code (ey, CPL_ERROR_NULL_INPUT);
525
526 /* Fill first guess */
527 cpl_array * parameters = cpl_array_new (7, CPL_TYPE_DOUBLE);
528 cpl_array_fill_window_invalid (parameters, 0, 7);
529 cpl_array_set (parameters, 2, 0);
530 cpl_array_set (parameters, 3, *x);
531 cpl_array_set (parameters, 4, *y);
532 cpl_array_set (parameters, 5, 3);
533 cpl_array_set (parameters, 6, 3);
534 CPLCHECK_MSG ("Error creating parameter table");
535
536 cpl_size x_med_min=(cpl_size)(*x)-size;
537 cpl_size x_med_max=(cpl_size)(*x)+size;
538 cpl_size y_med_min=(cpl_size)(*y)-size;
539 cpl_size y_med_max=(cpl_size)(*y)+size;
540 cpl_size nx = cpl_image_get_size_x (img);
541 cpl_size ny = cpl_image_get_size_y (img);
542 if (x_med_min <1) x_med_min=1;
543 if (y_med_min <1) y_med_min=1;
544 if (x_med_max > nx) x_med_max=nx;
545 if (y_med_max > ny) y_med_max=ny;
546
547 double med = cpl_image_get_median_window (img,
548 x_med_min, y_med_min,
549 x_med_max, y_med_max);
550 cpl_array_set (parameters, 0, med);
551 CPLCHECK_MSG ("Error getting median");
552
553 /* Fit Gaussian */
554 double rms=0.;
555
556 cpl_error_code fit_converged = cpl_fit_image_gaussian (img, NULL, (cpl_size)(*x), (cpl_size)(*y),
557 size, size, parameters,
558 NULL, NULL, &rms, NULL, NULL,
559 NULL, NULL, NULL, NULL);
560
561 if (fit_converged == CPL_ERROR_NONE)
562 {
563 /* Set back */
564 *x = cpl_array_get (parameters,3,NULL);
565 *y = cpl_array_get (parameters,4,NULL);
566 CPLCHECK_MSG ("Error reading fit result");
567
568 /* Check errors */
569 /* reject result if either: */
570 /* * peak is below 3*rms */
571 /* * sigma_x > 5 pixels */
572 /* * sigma_y > 5 pixels */
573 double A = cpl_array_get (parameters, 1, NULL);
574 double rho = cpl_array_get (parameters, 2, NULL);
575 double sx = cpl_array_get (parameters, 5, NULL);
576 double sy = cpl_array_get (parameters, 6, NULL);
577
578 // cf. Condon 1996, PASP 109:166
579 double cst = 2. * sqrt(2. * M_PI * (1. - rho*rho) * sx * sy) * rms / A;
580 *ex=cst*sx;
581 *ey=cst*sy;
582
583 if ( A < 0. ) {
584 // detection is just not significant
585 cpl_msg_info (cpl_func, "rejecting fit: x=%g, y=%g, SNR=%g, sx=%g, sy=%g",
586 *x, *y, A/(rms * 2.*M_PI*sx*sy*sqrt(1-rho*rho)), sx, sy);
587 *x = 0.;
588 *y = 0.;
589 *ex = -1.;
590 *ey = -1.;
591 }
592 CPLCHECK_MSG("Error checking significance of fit result");
593 } else {
594 cpl_msg_warning (cpl_func, "fit of pupil beacon did not converge");
595 *x = 0.;
596 *y = 0.;
597 *ex = -1.;
598 *ey = -1.;
599 cpl_error_reset();
600 }
601 CPLCHECK_MSG ("Pupil beacon fit failed");
602
603 /* Fill image with zero at the detected position */
604 //if (*x > 0. && *y > 0.) {
605 // cpl_image_set (img, (cpl_size)(*x), (cpl_size)(*y), 0.0);
606 //}
607 //CPLCHECK_MSG("Error setting peak to zero");
608
609 /* Delete */
610 FREE (cpl_array_delete, parameters);
611
613 return CPL_ERROR_NONE;
614}
615
616/*----------------------------------------------------------------------------*/
632/*----------------------------------------------------------------------------*/
633
634cpl_error_code gravi_acqcam_pupil_v2(cpl_image * mean_img, cpl_imagelist * acqcam_imglist, cpl_propertylist * header,
635 cpl_table * acqcam_table, cpl_propertylist * o_header,
636 gravi_data *static_param_data)
637{
639
640cpl_ensure_code(acqcam_imglist, CPL_ERROR_NULL_INPUT);
641cpl_ensure_code(header, CPL_ERROR_NULL_INPUT);
642cpl_ensure_code(acqcam_table, CPL_ERROR_NULL_INPUT);
643cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
644
645cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (acqcam_imglist, 0));
646cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (acqcam_imglist, 0));
647
648/*
649 * constant ury is the beginning of the upper part of the acquisition camera
650 * It delimitates where the pupil beacon images starts
651 * Not sure it works in full frame case (ny>1100)
652 */
653const cpl_size ury = (ny>1100) ? 1200 : 745;
654
655/*
656 * First step. The goal is to get the images pupil, cutted to keep only pupil frames (using ury).
657 * The output is the pupil, after kernel filtering with a Gaussian pdf the size of the pupil beacon
658 */
659 cpl_imagelist * pupilImage_filtered = cpl_imagelist_new ();
660 gravi_acqcam_clean_pupil_v2(acqcam_imglist, pupilImage_filtered, ury);
661 CPLCHECK_MSG("Cannot clean pupil");
662
663/*
664 * Second step. The goal is to select the frames with the blinking beacons on
665 * The same routines gets the frames where the beacons are off
666 * It uses theses frames to remove the background
667 */
668 cpl_array * good_frames= cpl_array_new(cpl_imagelist_get_size(pupilImage_filtered),CPL_TYPE_INT);
669 cpl_imagelist * pupilImage_onFrames = cpl_imagelist_new ();
670 gravi_acqcam_select_good_frames_v2(pupilImage_filtered,pupilImage_onFrames,good_frames);
671 CPLCHECK_MSG("Cannot find blinking pupil frames");
672
673/* Third step. The goal is to initialize the vectors which will be used to know where the diodes
674 * are on the detectors (hence the use of bivectors, which contains x and y positions
675 * the size of the arrays are defined using several fixed parameters:
676 * GRAVI_SPOT_NTEL = 4 corresponds to the 4 telescopes
677 * GRAVI_SPOT_NSPOT = 4 corresponds to the 4 diodes on each telescopes
678 * GRAVI_SPOT_NLENS = 4 corresponds to the 4 lenslets on the acquisition camera
679 * nrow_on= 122 (for example) which corresponds to the N frames where the diodes are ON */
680
681 cpl_size nrow = cpl_imagelist_get_size(pupilImage_filtered);
682 cpl_size nrow_on = cpl_imagelist_get_size(pupilImage_onFrames);
683
684 if ((nrow_on<1)||(nrow<4))
685 {
686 /* exiting with empty pupil table */
687 cpl_array_fill_window_int (good_frames, 0, nrow_on, 0);
688 gravi_acqcam_set_pupil_table_v2(acqcam_table, header, o_header, good_frames, NULL, NULL, NULL, NULL, NULL);
689 CPLCHECK_MSG("Cannot stor pupil offset values in OI_ACQ table");
690
691 cpl_imagelist_delete(pupilImage_onFrames);
692 cpl_imagelist_delete(pupilImage_filtered);
693 cpl_array_delete(good_frames);
694 cpl_msg_warning (cpl_func, "Cannot reduce the %lli frames (not enough frames)",nrow);
696 return CPL_ERROR_NONE;
697 }
698
699 cpl_vector * scale_vector = cpl_vector_new(GRAVI_SPOT_NTEL);
700 cpl_vector * focus_value = cpl_vector_new(GRAVI_SPOT_NTEL);
701 cpl_array * bad_frames_short= cpl_array_new(nrow_on,CPL_TYPE_INT);
702 cpl_array_fill_window_int (bad_frames_short, 0, nrow_on, 0);
703 cpl_bivector * diode_pos_subwindow = cpl_bivector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NTEL);
704 cpl_bivector ** diode_pos_telescope = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
705 cpl_bivector ** diode_pos_theoretical = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
706 cpl_bivector ** diode_pos_offset = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
707 cpl_imagelist ** pupilImage_shiftandadd = cpl_malloc (GRAVI_SPOT_NTEL * sizeof(cpl_imagelist *));
708
709 for (int n = 0 ; n < nrow_on; n++)
710 {
711 diode_pos_telescope[n] = cpl_bivector_new (GRAVI_SPOT_NSPOT*GRAVI_SPOT_NTEL);
712 diode_pos_theoretical[n] = cpl_bivector_new (GRAVI_SPOT_NFOCUS*GRAVI_SPOT_NSPOT*GRAVI_SPOT_NLENS*GRAVI_SPOT_NTEL);
713 diode_pos_offset[n] = cpl_bivector_new (GRAVI_SPOT_NTEL);
714 }
715
716 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
717 pupilImage_shiftandadd[tel] = cpl_imagelist_new();
718 CPLCHECK_MSG("Cannot initialized vectors");
719
720/*
721 * Fourth step. Here we get the position of the 4 subwindows (correspond to each lenslet), in pixel/detector space
722 */
723 gravi_acqcam_get_pup_ref_v2(header, diode_pos_subwindow);
724 CPLCHECK_MSG("Cannot get pupil reference values");
725
726/*
727 * Fifth step. Here we get the position of the 4 diodes of each 4 telescopes, in pixel/detector space
728 */
729 gravi_acqcam_get_diode_ref_v2(header, good_frames, scale_vector, diode_pos_telescope, nrow_on);
730 CPLCHECK_MSG("Cannot get telescope beacons position");
731
732/*
733 * Seventh step. Get a large array of diode vector position (diode_pos_theoretical)
734 * this for each telescope, diode, lenslet, and focal number.
735 */
736 gravi_acqcam_get_diode_theoretical_v2(diode_pos_subwindow ,diode_pos_telescope,
737 diode_pos_theoretical, nrow_on, ury);
738 CPLCHECK_MSG("Cannot get theoretical position of beacons");
739
740/*
741 * Eigth step. Test each focus position to see which one is best. Then, collapse all spot images according to that focus value.
742 * The result pupilImage_shiftandadd correspond to an image for each blinking frame and each telescope
743 */
744 gravi_acqcam_perform_shiftandadd_v2(pupilImage_onFrames, pupilImage_shiftandadd, good_frames, focus_value, diode_pos_theoretical,
745 diode_pos_offset, nrow_on);
746 CPLCHECK_MSG("Cannot find correct focus to shift and add pupil images");
747
748/*
749 * Nineth step. For each one of the pupilImage_shiftandadd, fit a gaussian to know where the pupil is
750 * The measurement is done as an offset with respect to the theroetical position of the beacons
751 */
752 gravi_acqcam_get_pupil_offset_v2(pupilImage_shiftandadd, bad_frames_short, diode_pos_offset , o_header, nrow_on);
753 CPLCHECK_MSG("Cannot get pupil offset values");
754
755/*
756 * Last step. Storing the value of the pupil position in a table
757 */
758 gravi_acqcam_set_pupil_table_v2(acqcam_table, header, o_header, good_frames, scale_vector, bad_frames_short, diode_pos_offset, focus_value, static_param_data);
759 CPLCHECK_MSG("Cannot stor pupil offset values in OI_ACQ table");
760
761 /*
762 * Bonus. Add collapsed image to mean image and imprint best position as a cross in image
763 */
764
765 cpl_image * pupilImage_onFrames_collapse=cpl_imagelist_collapse_create (pupilImage_onFrames);
766 cpl_image * pupilImage_onFrames_1frame=cpl_imagelist_get (pupilImage_onFrames,0);
767 gravi_image_replace_window (mean_img, pupilImage_onFrames_1frame, 1, ury, nx, ny, 1, 1);
768 gravi_acqcam_spot_imprint_v2(mean_img, diode_pos_offset, diode_pos_theoretical, ury);
769 CPLCHECK_MSG("Cannot modify mean_image");
770
771/* cleaning all arrays ( free memory ) */
772
773for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
774{
775 if (pupilImage_shiftandadd[tel] != NULL)
776 cpl_imagelist_delete( pupilImage_shiftandadd[tel]);
777}
778FREE (cpl_free, pupilImage_shiftandadd);
779
780
781 for (int n = 0 ; n < nrow_on; n++)
782 {
783 cpl_bivector_delete( diode_pos_telescope[n]);
784 cpl_bivector_delete( diode_pos_theoretical[n]);
785 cpl_bivector_delete( diode_pos_offset[n]);
786 }
787 FREE (cpl_free, diode_pos_telescope);
788 FREE (cpl_free, diode_pos_theoretical);
789 FREE (cpl_free, diode_pos_offset);
790
791 cpl_bivector_delete( diode_pos_subwindow);
792 cpl_image_delete(pupilImage_onFrames_collapse);
793 cpl_imagelist_delete(pupilImage_onFrames);
794 cpl_imagelist_delete(pupilImage_filtered);
795 cpl_array_delete(good_frames);
796 cpl_array_delete(bad_frames_short);
797 cpl_vector_delete(scale_vector);
798 cpl_vector_delete(focus_value);
799
800CPLCHECK_MSG("Failed at freeing memory");
802return CPL_ERROR_NONE;
803}
804
805/*----------------------------------------------------------------------------*/
816/*----------------------------------------------------------------------------*/
817
818cpl_error_code gravi_acqcam_clean_pupil_v2(cpl_imagelist * acqcam_imglist, cpl_imagelist * pupilImage_filtered, const cpl_size ury)
819{
821
822 cpl_ensure_code(acqcam_imglist, CPL_ERROR_NULL_INPUT);
823 cpl_ensure_code(pupilImage_filtered, CPL_ERROR_NULL_INPUT);
824 cpl_ensure_code(ury, CPL_ERROR_NULL_INPUT);
825
826 /* Number of row */
827 cpl_size nrow = cpl_imagelist_get_size(acqcam_imglist);
828
829 cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (acqcam_imglist, 0));
830 cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (acqcam_imglist, 0));
831
832 cpl_imagelist * pupilImage = cpl_imagelist_new ();
833
834 for (cpl_size n = 0; n < nrow ; n++)
835 {
836 cpl_image * small_img_tmp = cpl_imagelist_get (acqcam_imglist, n);
837 cpl_image * small_img = cpl_image_extract( small_img_tmp, 1, ury, nx, ny);
838 cpl_imagelist_set(pupilImage, small_img,n);
839 }
840
841 nx=cpl_image_get_size_x (cpl_imagelist_get (pupilImage, 0));
842 ny=cpl_image_get_size_y (cpl_imagelist_get (pupilImage, 0));
843
844 cpl_image* pupilImageBackground = cpl_imagelist_collapse_minmax_create (pupilImage,(cpl_size)(nrow/8),(cpl_size)(5*nrow/8));
845
846 cpl_imagelist_subtract_image (pupilImage,pupilImageBackground);
847 CPLCHECK_MSG("Failure to subtract background");
848
849 cpl_size Npref = 21;
850 cpl_matrix * kernel1 = cpl_matrix_new (Npref,Npref);
851 cpl_matrix * kernel2 = cpl_matrix_new (Npref,Npref);
852 for (cpl_size x =0;x< Npref; x++)
853 for (cpl_size y =0;y< Npref; y++)
854 {
855 double radius_square=(x-(Npref-1)/2)*(x-(Npref-1)/2)+(y-(Npref-1)/2)*(y-(Npref-1)/2);
856 cpl_matrix_set(kernel1,x,y,exp(- radius_square /30));
857 cpl_matrix_set(kernel2,x,y,exp(- radius_square /40));
858 }
859
860 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_mean (kernel1));
861 cpl_matrix_divide_scalar (kernel2, cpl_matrix_get_mean (kernel2));
862 cpl_matrix_subtract (kernel1,kernel2);
863 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_stdev (kernel1));
864 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_max (kernel1));
865
866 CPLCHECK_MSG("Error computing gaussian Kernel");
867
868 // zero-pad to match image and kernel size, plus extra buffer for spurious data at edges
869 cpl_size fft_nx = nx + Npref - 1;
870 cpl_size fft_ny = ny + Npref - 1;
871 cpl_size half_kernel_w = (Npref - 1) / 2;
872
873 // Copy kernel into image with zero-padding
874 cpl_matrix_resize(kernel1, 0, ny-1, 0, nx-1);
875 cpl_image *kernel_image = cpl_image_wrap_double(fft_nx, fft_ny, cpl_matrix_get_data(kernel1));
876
877 // Compute kernel FFT
878 cpl_image *kernel_fft = cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE_COMPLEX);
879 cpl_fft_image(kernel_fft, kernel_image, CPL_FFT_FORWARD);
880
881 // Create imagelist with zero-padded copies of the source images
882 cpl_imagelist *pupilImage_padded = cpl_imagelist_new();
883 for (cpl_size n = 0; n < nrow ; n++)
884 {
885 const cpl_image *small_img_tmp = cpl_imagelist_get_const(pupilImage, n);
886 cpl_image *padded_img_tmp = cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE);
887 cpl_image_copy(padded_img_tmp, small_img_tmp, 1, 1);
888 cpl_imagelist_set(pupilImage_padded, padded_img_tmp, n);
889 }
890
891 // Fourier-transform the images
892 cpl_imagelist *pupilImage_fft = cpl_imagelist_new();
893 for (cpl_size n = 0; n < nrow; n++)
894 cpl_imagelist_set(pupilImage_fft, cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE_COMPLEX), n);
895 cpl_fft_imagelist(pupilImage_fft, pupilImage_padded, CPL_FFT_FORWARD | CPL_FFT_FIND_MEASURE);
896
897 // Convolve with kernel
898 cpl_imagelist_multiply_image(pupilImage_fft, kernel_fft);
899
900 // Reverse FFT to get filtered images
901 cpl_imagelist *pupilImage_filtered_padded = cpl_imagelist_new();
902 for (cpl_size n = 0; n < nrow; n++)
903 cpl_imagelist_set(pupilImage_filtered_padded, cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE), n);
904 cpl_fft_imagelist(pupilImage_filtered_padded, pupilImage_fft, CPL_FFT_BACKWARD | CPL_FFT_FIND_MEASURE);
905
906 // Trim the padding from the image edges
907 for (cpl_size n = 0; n < nrow; n++) {
908 const cpl_image *pad_filt_img_tmp = cpl_imagelist_get_const(pupilImage_filtered_padded, n);
909 cpl_image *trimmed_img = cpl_image_extract(pad_filt_img_tmp, half_kernel_w+1, half_kernel_w+1, fft_nx-half_kernel_w, fft_ny-half_kernel_w);
910 cpl_imagelist_set(pupilImage_filtered, trimmed_img, n);
911 }
912
913 cpl_matrix_delete(kernel1);
914 cpl_matrix_delete(kernel2);
915
916 cpl_image_unwrap(kernel_image);
917 cpl_image_delete(kernel_fft);
918
919 cpl_image_delete(pupilImageBackground);
920 cpl_imagelist_delete(pupilImage);
921 cpl_imagelist_delete(pupilImage_fft);
922 cpl_imagelist_delete(pupilImage_padded);
923 cpl_imagelist_delete(pupilImage_filtered_padded);
924
925 CPLCHECK_MSG("Pupil Fitting V2 does not work");
927 return CPL_ERROR_NONE;
928}
929
930
931/*----------------------------------------------------------------------------*/
945/*----------------------------------------------------------------------------*/
946
947cpl_error_code gravi_acqcam_select_good_frames_v2(cpl_imagelist * pupilImage_filtered, cpl_imagelist * pupilImage_onFrames, cpl_array * good_frames)
948{
950 int nv =0;
951
952 cpl_ensure_code(pupilImage_filtered, CPL_ERROR_NULL_INPUT);
953 cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
954
955 /* Number of row , size of image */
956 cpl_size nrow = cpl_imagelist_get_size(pupilImage_filtered);
957 cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (pupilImage_filtered, 0));
958 cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (pupilImage_filtered, 0));
959
960 cpl_vector * pupil_max = cpl_vector_new(4);
961 cpl_vector * pupil_mediam = cpl_vector_new(nrow);
962 cpl_array * back_frames= cpl_array_new(nrow,CPL_TYPE_INT);
963
964 for (cpl_size n = 0; n < nrow; n++)
965 {
966
967 cpl_image * image = cpl_imagelist_get (pupilImage_filtered,n);
968 cpl_vector_set(pupil_max, 0, cpl_image_get_max_window (image, 1, 1, 250, ny));
969 cpl_vector_set(pupil_max, 1, cpl_image_get_max_window (image, 251, 1, 500, ny));
970 cpl_vector_set(pupil_max, 2, cpl_image_get_max_window (image, 501, 1, 750, ny));
971 cpl_vector_set(pupil_max, 3, cpl_image_get_max_window (image, 751, 1, nx, ny));
972 cpl_vector_set(pupil_mediam, n, cpl_vector_get_median (pupil_max));
973 /*
974 cpl_msg_info (cpl_func, "tototo image valeur %lli = %.2f",n, cpl_image_get(image, 101, 151, &nv));
975 cpl_msg_info (cpl_func, "tototo %.2f, %.2f, %.2f, %.2f", cpl_vector_get(pupil_max, 0),cpl_vector_get(pupil_max, 1),cpl_vector_get(pupil_max, 2),cpl_vector_get(pupil_max, 3));
976 cpl_msg_info (cpl_func, "Finding min brightness level for beacon ON %lli: %.5f ADU",n, cpl_vector_get(pupil_mediam, n));*/
977
978 }
979
980 cpl_vector * pupil_mediam_sort = cpl_vector_duplicate (pupil_mediam);
981 cpl_vector_sort (pupil_mediam_sort,CPL_SORT_ASCENDING);
982 cpl_vector * pupil_median_down = cpl_vector_extract (pupil_mediam_sort, 0, (nrow-1)/2, 1);
983 cpl_vector * pupil_median_up = cpl_vector_extract (pupil_mediam_sort, nrow/2, nrow-1, 1);
984 double threshold_up = cpl_vector_get_median(pupil_median_up);
985 double threshold_down = cpl_vector_get_median(pupil_median_down);
986 double threshold=(threshold_up+threshold_down)/2.0;
987
988 CPLCHECK_MSG("Could not get median maximum of beacon flux ");
989 cpl_msg_info (cpl_func, "Found threshold brightness level for beacon : %.5f ADU", threshold);
990
991 for (cpl_size n = 0; n < nrow; n++)
992 {
993 if (cpl_vector_get(pupil_mediam, n) >threshold *1.1)
994 cpl_array_set_int (good_frames,n,1);
995 else
996 cpl_array_set_int (good_frames,n,0);
997 if (cpl_vector_get(pupil_mediam, n) <threshold * 0.9)
998 cpl_array_set_int (back_frames,n,1);
999 else
1000 cpl_array_set_int (back_frames,n,0);
1001 }
1002
1003 cpl_size frames_background = 4;
1004 cpl_size n_goodFrames = 0;
1005
1006 CPLCHECK_MSG("Failed to find blincking pupil files");
1007
1008 for (cpl_size n = 0; n < nrow; n++)
1009 if (cpl_array_get_int(good_frames,n,&nv) == 1)
1010 {
1011 cpl_image * image = cpl_imagelist_get (pupilImage_filtered,n);
1012 /*cpl_msg_info (cpl_func, "Good frame image valeur %lli = %.2f",n, cpl_image_get(image, 101, 151, &nv));*/
1013
1014 cpl_size n_frames_background = 0;
1015 cpl_imagelist * pupilImage_offFrames = cpl_imagelist_new ();
1016 for (cpl_size b = 1; b < nrow; b++)
1017 {
1018 if (n+b < nrow)
1019 if ((cpl_array_get_int(back_frames,n+b,&nv) == 1)&(n_frames_background<frames_background))
1020 {
1021 cpl_image * image_background = cpl_imagelist_get (pupilImage_filtered,n+b);
1022 cpl_image * image_background_copy = cpl_image_duplicate (image_background);
1023 cpl_imagelist_set(pupilImage_offFrames, image_background_copy, n_frames_background);
1024 n_frames_background = n_frames_background + 1;
1025 }
1026 if (n-b >= 0)
1027 if ((cpl_array_get_int(back_frames,n-b,&nv) == 1)&(n_frames_background<frames_background))
1028 {
1029 cpl_image * image_background = cpl_imagelist_get (pupilImage_filtered,n-b);
1030 cpl_image * image_background_copy = cpl_image_duplicate (image_background);
1031 cpl_imagelist_set(pupilImage_offFrames, image_background_copy, n_frames_background);
1032 n_frames_background = n_frames_background + 1;
1033 }
1034 }
1035
1036 cpl_image * image_background_mean = cpl_imagelist_collapse_create (pupilImage_offFrames);
1037 cpl_image * image_background_subtracted = cpl_image_subtract_create (image, image_background_mean);
1038 cpl_imagelist_set(pupilImage_onFrames, image_background_subtracted, n_goodFrames);
1039
1040 cpl_image_delete(image_background_mean);
1041 cpl_imagelist_delete(pupilImage_offFrames);
1042
1043 n_goodFrames ++;
1044 }
1045 cpl_msg_info(cpl_func, "Found %lli frames with pupil beacon ON (over %lli frames)",n_goodFrames, nrow);
1046
1047 cpl_vector_delete(pupil_max);
1048 cpl_vector_delete(pupil_mediam);
1049 cpl_vector_delete(pupil_median_up);
1050 cpl_vector_delete(pupil_median_down);
1051 cpl_vector_delete(pupil_mediam_sort);
1052 cpl_array_delete(back_frames);
1053
1054 CPLCHECK_MSG("Pupil Selection of good files failed");
1056 return CPL_ERROR_NONE;
1057}
1058
1059/*----------------------------------------------------------------------------*/
1071/*----------------------------------------------------------------------------*/
1072
1073cpl_error_code gravi_acqcam_get_pup_ref_v2 (cpl_propertylist * header, cpl_bivector * diode_pos_subwindow)
1074{
1076 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1077 cpl_ensure_code (diode_pos_subwindow, CPL_ERROR_NULL_INPUT);
1078
1079 cpl_size nsx = 0, nsy = 0, sx = 0, sy = 0;
1080
1081 cpl_vector * x_pos_subwindow = cpl_bivector_get_x (diode_pos_subwindow);
1082 cpl_vector * y_pos_subwindow = cpl_bivector_get_y (diode_pos_subwindow);
1083
1084
1085 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1086 {
1087 /* If sub-windowing, we read the sub-window size and
1088 * the sub-window start for pupil */
1089 if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
1090 char name[90];
1091
1092 nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
1093 sprintf (name, "ESO DET1 FRAM%d STRX", 3*GRAVI_SPOT_NTEL + tel + 1);
1094 sx = cpl_propertylist_get_int (header, name);
1095
1096 nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY");
1097 sprintf (name, "ESO DET1 FRAM%d STRY", 3*GRAVI_SPOT_NTEL + tel + 1);
1098 sy = cpl_propertylist_get_int (header, name);
1099
1100 CPLCHECK_MSG ("Cannot get sub-windowing parameters");
1101 }
1102
1103 cpl_msg_info (cpl_func,"sub-window pupil %i sx= %lld sy = %lld", tel, sx, sy);
1104
1105 /* Read the sub-apperture reference positions
1106 * Converted to accound for sub-windowing
1107 * In vector convention (start at 0,0) */
1108
1109 for (int lens = 0; lens < GRAVI_SPOT_NLENS ; lens++) {
1110 double xv = gravi_pfits_get_ptfc_acqcam (header, lens*GRAVI_SPOT_NTEL + tel + 1);
1111 double yv = gravi_pfits_get_ptfc_acqcam (header, lens*GRAVI_SPOT_NTEL + tel + 17);
1112 cpl_vector_set( x_pos_subwindow, lens*GRAVI_SPOT_NTEL+tel , xv - (sx - tel*nsx) - 1);
1113 cpl_vector_set( y_pos_subwindow, lens*GRAVI_SPOT_NTEL+tel , yv - (sy - 3*nsy) - 1);
1114 /* cpl_msg_info (cpl_func,"pupil %lli subC %i = %10.4f,%10.4f",
1115 tel, lens, x_diode_subwindow[lens][tel], y_diode_subwindow[lens][tel]);*/
1116 CPLCHECK_MSG ("Cannot get pupil reference position");
1117 }
1118 }
1119
1121 return CPL_ERROR_NONE;
1122}
1123
1124
1125/*----------------------------------------------------------------------------*/
1140/*----------------------------------------------------------------------------*/
1141
1142
1143cpl_error_code gravi_acqcam_get_diode_ref_v2 (cpl_propertylist * header,
1144 cpl_array * good_frames,
1145 cpl_vector * scale_vector,
1146 cpl_bivector ** diode_pos_telescope,
1147 int nrow_on)
1148{
1150 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1151 cpl_ensure_code (good_frames, CPL_ERROR_NULL_INPUT);
1152
1153 int nv = 0;
1154 double dx,dy,scale;
1155 double parang1 = cpl_propertylist_get_double(header, "ESO ISS PARANG START");
1156 double parang2 = cpl_propertylist_get_double(header, "ESO ISS PARANG END");
1157 CPLCHECK_MSG("Cannot determine parallactic angle");
1158
1159 /* check the feed mode */
1160 int sts_mode = 0;
1161 if ( !strcmp(gravi_pfits_get_feed (header), "DUAL_STS") ) {
1162 sts_mode = DUAL_STS;
1163 }
1164 else{
1165 sts_mode = SINGLE_STS;
1166 }
1167
1168 cpl_size nrow = cpl_array_get_size(good_frames);
1169 if (fabs(cpl_array_get_mean(good_frames)*nrow -nrow_on) > 1e-2)
1170 cpl_msg_error (cpl_func, "Ratio of blinking frames different %f from %i",cpl_array_get_mean(good_frames)*nrow,nrow_on);
1171
1172 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1173 {
1174
1175 /* Get the telescope name and ID */
1176 const char * telname = gravi_conf_get_telname (tel, header);
1177
1178 /* Check telescope name */
1179 if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
1180 cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
1181 CPLCHECK_MSG("Cannot get telescope names");
1182
1183 /* Hardcoded theoretical positions in mm */
1184
1185 /* If UTs or ATs, select scaling, rotation, and spacing */
1186 if (telname[0] == 'U') {
1187 // FE 2019-05-15 replaced with median measured for the whole 2017-2018
1188 // Galactic Center data set, which should be the best information at hand
1189 // cpl_vector_set (output, GRAVI_SPOT_SCALE, 16.225);
1190 if ( tel == 0 ) scale = 16.83;
1191 if ( tel == 1 ) scale = 17.42;
1192 if ( tel == 2 ) scale = 16.85;
1193 if ( tel == 3 ) scale = 17.41;
1194 // below information could be calculated from diode position in header
1195 // this would also give the 45 offset angle introduced above to calculate
1196 // the spot angle
1197 dx = 0.363;
1198 dy = 0.823;
1199 } else if (telname[0] == 'A') {
1200 // FE 2019-05-15 maybe we should also update the AT numbers
1201 scale = 73.0154;
1202 // EW, FE 2019-09-11: short and long side of AT beacons are
1203 // flipped when compared to UTs
1204 dx = 0.158;
1205 dy = 0.122;
1206 } else
1207 return cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1208 "Cannot get telescope name");
1209 CPLCHECK_MSG("Cannot get telescope scale");
1210 cpl_vector_set(scale_vector,tel,scale);
1211
1212 /* We see that fitting the spot angle fails if one spot is fainter than reflection from brightest spot
1213 resulting in wrong pupil_x,y, therefore we calculate angle from header information
1214 by design of telescope spiders. This model angle should be
1215 angle = north_angle + paralactic angle + 45
1216 Remark: checking the angles for the GC easter flare night, we get an offset of 45.75 degrees */
1217 double northangle = 0;
1218 double zenithangle = 0;
1219 if (sts_mode == SINGLE_STS) {
1220 northangle = gravi_pfits_get_northangle_acqcam(header, tel);
1221 CPLCHECK_MSG("Cannot determine field angle");
1222 }
1223
1224
1225 cpl_size n_on= 0;
1226 for (cpl_size n =0; n < nrow; n++)
1227 {
1228 /*cpl_msg_info (cpl_func, "Totos: %lli, %.2f",cpl_array_get_int(good_frames,n,&nv),cpl_array_get_mean(good_frames));*/
1229 if (cpl_array_get_int(good_frames,n,&nv) ==1)
1230 {
1231 double padif = parang2 - parang1;
1232 if (padif > 180)
1233 padif -= 360;
1234 if (padif < -180)
1235 padif += 360;
1236 double parang= parang1 + padif*n/(nrow-1);
1237
1238 // TODO: FE handle case of calibration unit //
1239 double angle = 0;
1240 if (sts_mode == DUAL_STS) {
1241 /* northangle = gravi_pfits_get_northangle_beamb_acqcam(header, tel, n, nrow); */
1242 zenithangle = gravi_pfits_get_zenithangle_beamb_acqcam(header, tel, n, nrow);
1243 angle = zenithangle - 45.;
1244 } else {
1245 angle = northangle + parang + 45.;
1246 }
1247
1248 if (angle < 0) angle += 180;
1249 if (angle > 180) angle -= 180;
1250 double cang = cos(angle * CPL_MATH_RAD_DEG) * scale;
1251 double sang = sin(angle * CPL_MATH_RAD_DEG) * scale;
1252
1253 /* Diode arrangement */
1254
1255 cpl_vector * x_pos_telescope = cpl_bivector_get_x (diode_pos_telescope[n_on]);
1256 cpl_vector * y_pos_telescope = cpl_bivector_get_y (diode_pos_telescope[n_on]);
1257
1258 cpl_vector_set(x_pos_telescope, 0*GRAVI_SPOT_NTEL+tel, -cang * dx + sang * dy);
1259 cpl_vector_set(x_pos_telescope, 1*GRAVI_SPOT_NTEL+tel, -cang * dx - sang * dy);
1260 cpl_vector_set(x_pos_telescope, 2*GRAVI_SPOT_NTEL+tel, cang * dx - sang * dy);
1261 cpl_vector_set(x_pos_telescope, 3*GRAVI_SPOT_NTEL+tel, cang * dx + sang * dy);
1262
1263 cpl_vector_set(y_pos_telescope, 0*GRAVI_SPOT_NTEL+tel, cang * dy + sang * dx);
1264 cpl_vector_set(y_pos_telescope, 1*GRAVI_SPOT_NTEL+tel, -cang * dy + sang * dx);
1265 cpl_vector_set(y_pos_telescope, 2*GRAVI_SPOT_NTEL+tel, -cang * dy - sang * dx);
1266 cpl_vector_set(y_pos_telescope, 3*GRAVI_SPOT_NTEL+tel, cang * dy - sang * dx);
1267
1268 n_on ++;
1269
1270 if (n_on == 10)
1271 {
1272 cpl_msg_info (cpl_func, "angle, parang %.2f, %.2f",angle,parang);
1273 cpl_msg_info (cpl_func, "dx and dy %.2f, %.2f",cang,sang);
1274 }
1275
1276
1277 }
1278 }
1279
1280
1281 CPLCHECK_MSG("Cannot determine diode position in telescope space");
1282 }
1283
1285 return CPL_ERROR_NONE;
1286}
1287
1288/*----------------------------------------------------------------------------*/
1306/*----------------------------------------------------------------------------*/
1307
1308cpl_error_code gravi_acqcam_get_diode_theoretical_v2(cpl_bivector * diode_pos_subwindow,
1309 cpl_bivector ** diode_pos_telescope,
1310 cpl_bivector ** diode_pos_theoretical,
1311 cpl_size nrow_on, int ury)
1312{
1313
1315 cpl_ensure_code (diode_pos_subwindow, CPL_ERROR_NULL_INPUT);
1316 cpl_ensure_code (diode_pos_telescope, CPL_ERROR_NULL_INPUT);
1317 cpl_ensure_code (diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1318
1319 cpl_vector * x_pos_subwindow = cpl_bivector_get_x (diode_pos_subwindow);
1320 cpl_vector * y_pos_subwindow = cpl_bivector_get_y (diode_pos_subwindow);
1321
1322 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1323 {
1324 double x_lenslet_mean=0.0;
1325 double y_lenslet_mean=0.0;
1326 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1327 {
1328 x_lenslet_mean+=cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)/GRAVI_SPOT_NLENS;
1329 y_lenslet_mean+=cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)/GRAVI_SPOT_NLENS;
1330 }
1331
1332 cpl_msg_info (cpl_func, "Reference pixel position for tel %d : X = %.2f, Y= %.2f", tel, x_lenslet_mean, y_lenslet_mean);
1333
1334 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1335 for (int focus = 0 ; focus < GRAVI_SPOT_NFOCUS; focus++)
1336 {
1337 double xlenslet=cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)+
1338 (cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)-x_lenslet_mean)*gravi_acqcam_defocus_scaling(focus);
1339 double ylenslet=cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)+
1340 (cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)-y_lenslet_mean)*gravi_acqcam_defocus_scaling(focus);
1341
1342 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1343 for (cpl_size n_on = 0 ; n_on < nrow_on; n_on++)
1344 {
1345
1346 cpl_vector * x_pos_telescope = cpl_bivector_get_x (diode_pos_telescope[n_on]);
1347 cpl_vector * y_pos_telescope = cpl_bivector_get_y (diode_pos_telescope[n_on]);
1348 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n_on]);
1349 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n_on]);
1350
1351 double x_theo = xlenslet + cpl_vector_get(x_pos_telescope,spot*GRAVI_SPOT_NTEL+tel);
1352 double y_theo = ylenslet + cpl_vector_get(y_pos_telescope,spot*GRAVI_SPOT_NTEL+tel)- ury +1;
1353
1354 cpl_vector_set(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel,x_theo);
1355 cpl_vector_set(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel,y_theo);
1356 }
1357 }
1358 }
1359
1360CPLCHECK_MSG("Cannot compute the theoretical coordinates of spots");
1362return CPL_ERROR_NONE;
1363
1364}
1365
1366/*----------------------------------------------------------------------------*/
1376/*----------------------------------------------------------------------------*/
1377
1378
1379cpl_error_code gravi_acqcam_spot_imprint_v2(cpl_image *mean_img, cpl_bivector ** diode_pos_offset, cpl_bivector ** diode_pos_theoretical, int ury)
1380{
1382
1383 cpl_ensure_code(mean_img, CPL_ERROR_NULL_INPUT);
1384 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1385 cpl_ensure_code(diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1386
1387 cpl_size nx = cpl_image_get_size_x (mean_img);
1388 cpl_size ny = cpl_image_get_size_y (mean_img);
1389
1390 /* imprint position for first image */
1391 int n_on = 0;
1392 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n_on]);
1393 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n_on]);
1394 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n_on]);
1395 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n_on]);
1396
1397 /* imprint position with no defocus */
1398 int focus = GRAVI_SPOT_NFOCUS/2;
1399
1400 /* Loop on diode and appertures */
1401 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1402 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1403 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1404 {
1405 double x_off = cpl_vector_get(x_pos_offset,tel);
1406 double y_off = cpl_vector_get(y_pos_offset,tel)+ury;
1407 double x_diode = cpl_vector_get(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1408 double y_diode = cpl_vector_get(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1409 cpl_size xf = roundl(x_diode + x_off) + 1;
1410 cpl_size yf = roundl(y_diode + y_off) + 1;
1411 if (xf < 2 || xf > nx-2 || yf < 2 || yf > ny-2) continue;
1412 cpl_image_set (mean_img, xf, yf, 0);
1413 cpl_image_set (mean_img, xf-1, yf, 0);
1414 cpl_image_set (mean_img, xf+1, yf, 0);
1415 cpl_image_set (mean_img, xf, yf+1, 0);
1416 cpl_image_set (mean_img, xf, yf-1, 0);
1417 CPLCHECK ("Cannot imprint cross in image");
1418 }
1419
1420 CPLCHECK ("Cannot imprint cross in image");
1422 return CPL_ERROR_NONE;
1423}
1424
1425
1426/*----------------------------------------------------------------------------*/
1457/*----------------------------------------------------------------------------*/
1458
1459
1460cpl_error_code gravi_acqcam_perform_shiftandadd_v2(cpl_imagelist * pupilImage_onFrames, cpl_imagelist ** pupilImage_shiftandadd, cpl_array * good_frames,
1461 cpl_vector * focus_value,
1462 cpl_bivector ** diode_pos_theoretical ,
1463 cpl_bivector ** diode_pos_offset, cpl_size nrow_on)
1464{
1466
1467 cpl_ensure_code(pupilImage_onFrames, CPL_ERROR_NULL_INPUT);
1468 cpl_ensure_code(pupilImage_shiftandadd, CPL_ERROR_NULL_INPUT);
1469 cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
1470 cpl_ensure_code(focus_value, CPL_ERROR_NULL_INPUT);
1471 cpl_ensure_code(diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1472 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1473 cpl_ensure_code(nrow_on, CPL_ERROR_NULL_INPUT);
1474
1475 if (cpl_imagelist_get_size(pupilImage_onFrames) != nrow_on)
1476 cpl_msg_error (cpl_func, "Problem with number of blinking frames");
1477
1478 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1479 {
1480 cpl_vector * focus_max = cpl_vector_new(GRAVI_SPOT_NFOCUS);
1481 for (cpl_size focus = 0 ; focus < GRAVI_SPOT_NFOCUS; focus++)
1482 {
1483 cpl_imagelist * pupilImage_shifted= cpl_imagelist_new();
1484 cpl_size n_images = 0;
1485 for (cpl_size n = 0 ; n < nrow_on; n++)
1486 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1487 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1488 {
1489 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n]);
1490 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n]);
1491
1492 double x_theo=cpl_vector_get(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1493 double y_theo=cpl_vector_get(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1494
1495 cpl_size x=(int) (x_theo + 0.5);
1496 cpl_size y=(int) (y_theo+ 0.5);
1497 cpl_image * verysmall_img_tmp = cpl_imagelist_get (pupilImage_onFrames, n);
1498 cpl_image * verysmall_img = gravi_image_extract( verysmall_img_tmp, x+1-GRAVI_SPOT_NSEARCH, y+1-GRAVI_SPOT_NSEARCH, x+1+GRAVI_SPOT_NSEARCH, y+1+GRAVI_SPOT_NSEARCH);
1499 cpl_imagelist_set(pupilImage_shifted, verysmall_img, n_images);
1500 n_images ++;
1501 }
1502
1503 cpl_image * image_mean = cpl_imagelist_collapse_create (pupilImage_shifted);
1504 cpl_vector_set(focus_max, focus, cpl_image_get_max(image_mean));
1505 cpl_imagelist_delete(pupilImage_shifted);
1506 cpl_image_delete(image_mean);
1507 }
1508
1509 int focus_max_pos = cpl_vector_get_maxpos (focus_max);
1510 cpl_msg_info (cpl_func, "Focus value for telescope %d : F = %.2f pixels", tel, GRAVI_SPOT_SWINDOW * gravi_acqcam_defocus_scaling(focus_max_pos));
1511 cpl_vector_set(focus_value,tel,gravi_acqcam_defocus_scaling(focus_max_pos)*100.);
1512 CPLCHECK_MSG("Cannot find optimum focus position");
1513
1514 for (cpl_size n = 0 ; n < nrow_on; n++)
1515 {
1516 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n]);
1517 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n]);
1518 cpl_imagelist * pupilImage_shifted= cpl_imagelist_new();
1519 cpl_vector * offs_x = cpl_vector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NSPOT);
1520 cpl_vector * offs_y = cpl_vector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NSPOT);
1521 cpl_size n_images = 0;
1522 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1523 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1524 {
1525 double x_theo=cpl_vector_get(x_pos_theoretical,((focus_max_pos*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1526 double y_theo=cpl_vector_get(y_pos_theoretical,((focus_max_pos*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1527 int x= (int) (x_theo+ 0.5);
1528 int y= (int) (y_theo + 0.5);
1529 cpl_image * verysmall_img_tmp = cpl_imagelist_get (pupilImage_onFrames, n);
1530 cpl_image * verysmall_img = gravi_image_extract( verysmall_img_tmp, x+1-GRAVI_SPOT_NSEARCH, y+1-GRAVI_SPOT_NSEARCH, x+1+GRAVI_SPOT_NSEARCH, y+1+GRAVI_SPOT_NSEARCH);
1531 /*if (n==10)
1532 {
1533 cpl_msg_info (cpl_func, "Cutting at position : X/Y = %lli/%lli, [%lli,%lli,%lli]", x,y,lens,spot,tel);
1534 cpl_msg_info (cpl_func, "image_tmp %lli shifted : X/Y = %lli/%lli/%lli 15/10 %.2f", n,tel,lens,spot,cpl_image_get(verysmall_img_tmp, 11, 16, &nv));
1535 cpl_msg_info (cpl_func, "image %lli shifted : X/Y = %lli/%lli/%lli 15/10 %.2f", nrow_on,tel,lens,spot,cpl_image_get(verysmall_img, 11, 16, &nv));
1536 }*/
1537 /*if (n == 10)
1538 cpl_msg_info (cpl_func, "image %lli shifted : X/Y = %lli/%lli/%lli 15/10 %.2f -- off %.2f", n,tel,lens,spot,cpl_image_get(verysmall_img, 11, 16, &nv),x_diode_theo_fine[n][focus_max_pos][spot][lens][tel]);*/
1539 cpl_imagelist_set(pupilImage_shifted, verysmall_img, n_images);
1540 cpl_vector_set (offs_x, n_images, (x-x_theo));
1541 cpl_vector_set (offs_y, n_images, (y-y_theo));
1542
1543 n_images ++;
1544 }
1545
1546 cpl_bivector * offs = cpl_bivector_wrap_vectors (offs_x, offs_y);
1547
1548 double ppos_x, ppos_y;
1549 cpl_image** cpl_image_combined = cpl_geom_img_offset_saa (pupilImage_shifted,offs,CPL_KERNEL_DEFAULT,
1550 0,0,CPL_GEOM_INTERSECT,&ppos_x,&ppos_y);
1551 CPLCHECK_MSG("Cannot do fine shift and add");
1552
1553 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n]);
1554 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n]);
1555
1556 cpl_image * image_mean = cpl_image_extract (cpl_image_combined[0], 4, 4, GRAVI_SPOT_NSEARCH*2-3, GRAVI_SPOT_NSEARCH*2-3);
1557
1558 cpl_vector_set(x_pos_offset, tel, ppos_x+3);
1559 cpl_vector_set(y_pos_offset, tel, ppos_y+3);
1560
1561 cpl_imagelist_set(pupilImage_shiftandadd[tel], image_mean, n);
1562 CPLCHECK_MSG("Cannot add shift and added image to imagelist");
1563
1564 if (cpl_image_combined[0] != NULL) cpl_image_delete(cpl_image_combined[0]);
1565 if (cpl_image_combined[1] != NULL) cpl_image_delete(cpl_image_combined[1]);
1566 cpl_free(cpl_image_combined);
1567
1568 cpl_bivector_delete(offs);
1569 cpl_imagelist_delete(pupilImage_shifted);
1570 }
1571 cpl_vector_delete(focus_max);
1572 }
1573
1574 CPLCHECK_MSG("Fail at running the shift and add ESO algorithm");
1576
1577 return CPL_ERROR_NONE;
1578}
1579
1580
1581/*----------------------------------------------------------------------------*/
1599/*----------------------------------------------------------------------------*/
1600
1601
1602cpl_error_code gravi_acqcam_get_pupil_offset_v2(cpl_imagelist ** pupilImage_shiftandadd,
1603 cpl_array * bad_frames_short,
1604 cpl_bivector ** diode_pos_offset, cpl_propertylist * o_header,
1605 cpl_size nrow_on)
1606{
1608
1609 cpl_ensure_code(pupilImage_shiftandadd, CPL_ERROR_NULL_INPUT);
1610 cpl_ensure_code(bad_frames_short, CPL_ERROR_NULL_INPUT);
1611 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1612 cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
1613 cpl_ensure_code(nrow_on, CPL_ERROR_NULL_INPUT);
1614
1615 int nv =0;
1616
1617 if (cpl_array_get_size(bad_frames_short)!= nrow_on)
1618 cpl_msg_error (cpl_func, "Problem with number of blinking frames");
1619
1620 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1621 {
1622 double previous_xsc=0.0;
1623 double previous_ysc=0.0;
1624 double flux_sum=0.0;
1625 int n_sum=0;
1626
1627 for (cpl_size n = 0 ; n < nrow_on; n++)
1628 {
1629 cpl_image * image = cpl_imagelist_get (pupilImage_shiftandadd[tel], n);
1630 CPLCHECK_MSG("Failure at getting saa image");
1631
1632 cpl_size px, py;
1633 cpl_image_get_maxpos (image, &px, &py);
1634 CPLCHECK_MSG("Failure to read maximum position");
1635
1636 double flux_max= cpl_image_get(image,px,py,&nv);
1637 double std = cpl_image_get_stdev (image);
1638
1639
1640 double xsc = (double) px;
1641 double ysc = (double) py;
1642 double exsc=3., eysc=3.;
1643 cpl_size size = 15;
1644
1645 gravi_acq_fit_gaussian (image, &xsc, &ysc, &exsc, &eysc, size);
1646 CPLCHECK_MSG("Failure at fitting pupil spot");
1647
1648 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n]);
1649 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n]);
1650 CPLCHECK_MSG("Failure at reading bivector offsets");
1651
1652 double x_offset_final= cpl_vector_get(x_pos_offset, tel) + xsc-GRAVI_SPOT_NSEARCH-1;
1653 double y_offset_final= cpl_vector_get(y_pos_offset, tel) + ysc-GRAVI_SPOT_NSEARCH-1;
1654
1655 cpl_vector_set(x_pos_offset, tel, x_offset_final);
1656 cpl_vector_set(y_pos_offset, tel, y_offset_final);
1657
1658 int previous;
1659 if ((flux_max/std < 4)||(exsc<0))
1660 {
1661 previous = cpl_array_get_int (bad_frames_short,n,&nv);
1662 cpl_array_set_int (bad_frames_short,n,previous|(1<<tel));
1663 } else {
1664 if (n>0)
1665 {
1666 double distance=(xsc-previous_xsc)*(xsc-previous_xsc)+(ysc-previous_ysc)*(ysc-previous_ysc);
1667 if (distance > 10*10)
1668 {
1669 previous = cpl_array_get_int (bad_frames_short,n-1,&nv);
1670 cpl_array_set_int (bad_frames_short,n-1,previous|(1<<(tel+4)));
1671 previous = cpl_array_get_int (bad_frames_short,n,&nv);
1672 cpl_array_set_int (bad_frames_short,n,previous|(1<<(tel+4)));
1673 }
1674 }
1675 flux_sum+=flux_max;
1676 n_sum+=1;
1677 }
1678 previous_xsc=xsc;
1679 previous_ysc=ysc;
1680 CPLCHECK_MSG("Failure at testing bad pupil measurments");
1681 }
1682
1683 char qc_name[100];
1684 double flux_mean=flux_sum/(1e-7+n_sum);
1685 /* add QC parameters with pupil beacon flux */
1686 sprintf(qc_name, "ESO QC ACQ PUP%i SPOT_FLUX", tel + 1);
1687 cpl_msg_info(cpl_func, "%s = %f", qc_name,flux_mean);
1688 cpl_propertylist_update_double(o_header, qc_name, flux_mean);
1689 cpl_propertylist_set_comment(o_header, qc_name,"[ADU] mean spot flux");
1690
1691 CPLCHECK_MSG("Failing to store QC value");
1692 }
1693
1694 CPLCHECK_MSG("Failure at fitting gaussian on diode position");
1695
1697 return CPL_ERROR_NONE;
1698}
1699
1700
1701/*----------------------------------------------------------------------------*/
1723/*----------------------------------------------------------------------------*/
1724
1725
1726
1727cpl_error_code gravi_acqcam_set_pupil_table_v2(cpl_table * acqcam_table, cpl_propertylist * header, cpl_propertylist * o_header, cpl_array * good_frames,
1728 cpl_vector* scale_vector, cpl_array * bad_frames_short, cpl_bivector ** diode_pos_offset, cpl_vector * focus_value, gravi_data *static_param_data)
1729{
1730
1732int nv =0;
1733
1734cpl_ensure_code(acqcam_table, CPL_ERROR_NULL_INPUT);
1735cpl_ensure_code(header, CPL_ERROR_NULL_INPUT);
1736cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
1737cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
1738cpl_size nrow = cpl_array_get_size(good_frames);
1739
1740/* Pupil positions array */
1741gravi_table_new_column(acqcam_table, "PUPIL_NSPOT", NULL, CPL_TYPE_INT);
1742gravi_table_new_column(acqcam_table, "PUPIL_X", "pix", CPL_TYPE_DOUBLE);
1743gravi_table_new_column(acqcam_table, "PUPIL_Y", "pix", CPL_TYPE_DOUBLE);
1744gravi_table_new_column(acqcam_table, "PUPIL_Z", "pix", CPL_TYPE_DOUBLE);
1745gravi_table_new_column(acqcam_table, "PUPIL_R", "deg", CPL_TYPE_DOUBLE);
1746gravi_table_new_column(acqcam_table, "PUPIL_U", "m", CPL_TYPE_DOUBLE);
1747gravi_table_new_column(acqcam_table, "PUPIL_V", "m", CPL_TYPE_DOUBLE);
1748gravi_table_new_column(acqcam_table, "PUPIL_W", "m", CPL_TYPE_DOUBLE);
1749gravi_table_new_column(acqcam_table, "OPD_PUPIL", "m", CPL_TYPE_DOUBLE);
1750
1751double sobj_x = gravi_pfits_get_sobj_x(header);
1752double sobj_y = gravi_pfits_get_sobj_y(header);
1753CPLCHECK_MSG("Cannot determine SOBJ X and Y");
1754
1755if (cpl_array_get_mean(good_frames)*nrow<.5)
1756{
1757 /* no good frames */
1758 /* filling the table with zeros */
1759 cpl_msg_info(cpl_func,"No good pupil images ==> no good pupil position available");
1760
1761 for (int tel = 0; tel < GRAVI_SPOT_NTEL; tel++)
1762 for (cpl_size row = 0; row < nrow; row++)
1763 {
1764 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1765 CPLCHECK_MSG("Cannot put 0 in NSPOT data in the ACQ PUPIL");
1766 }
1767
1769 return CPL_ERROR_NONE;
1770}
1771
1772 /* looks like the table will not be empty. Checking additional parameters*/
1773cpl_ensure_code(scale_vector, CPL_ERROR_NULL_INPUT);
1774cpl_ensure_code(bad_frames_short, CPL_ERROR_NULL_INPUT);
1775cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1776cpl_ensure_code(focus_value, CPL_ERROR_NULL_INPUT);
1777cpl_ensure_code(static_param_data, CPL_ERROR_NULL_INPUT);
1778
1779cpl_size nrow_on = cpl_array_get_size(bad_frames_short);
1780
1781if (fabs(cpl_array_get_mean(good_frames)*nrow -nrow_on) > 1e-2)
1782 cpl_msg_error (cpl_func, "Ratio of blinking frames different %f from %lli",cpl_array_get_mean(good_frames)*nrow,nrow_on);
1783
1784
1785
1786for (int tel = 0; tel < GRAVI_SPOT_NTEL; tel++)
1787{
1788 int n_on=0;
1789
1790 /* Get the conversion angle xy to uv in [rad] */
1791 double northangle = gravi_pfits_get_northangle_acqcam(header, tel);
1792 double cfangle = cos(northangle * CPL_MATH_RAD_DEG);
1793 double sfangle = sin(northangle * CPL_MATH_RAD_DEG);
1794 double scale = cpl_vector_get(scale_vector,tel);
1795 CPLCHECK_MSG("Cannot determine field angle");
1796
1797 /* values which are valids for all frames */
1798 double r_shift = 0.0;
1799 double z_shift = cpl_vector_get(focus_value,tel)/GRAVI_SPOT_SWINDOW;
1800 double w_shift = gravi_acqcam_z2meter(z_shift, static_param_data);
1801
1802 /*values that does change */
1803 double x_shift = 0.0;
1804 double y_shift = 0.0;
1805 double u_shift = 0.0;
1806 double v_shift = 0.0;
1807 double opd_pupil = 0.0;
1808
1809 /* buffers to average values */
1810 double x_shift_sum = 0.0;
1811 double y_shift_sum = 0.0;
1812 double u_shift_sum = 0.0;
1813 double v_shift_sum = 0.0;
1814 int n_sum = 0;
1815 int n_bad_snr = 0;
1816 int n_bad_distance = 0;
1817
1818 for (cpl_size row = 0; row < nrow; row++)
1819 {
1820 /*cpl_msg_info (cpl_func, "Row number %lli n %lli",row,n_on);*/
1821
1822 if (cpl_array_get_int(good_frames,row,&nv) != 1)
1823 {
1824 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1825 CPLCHECK_MSG("Cannot put 0 in NSPOT data in the ACQ PUPIL");
1826 }
1827 else
1828 {
1829
1830 int is_it_bad=cpl_array_get_int (bad_frames_short,n_on,&nv);
1831
1832 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n_on]);
1833 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n_on]);
1834 CPLCHECK_MSG("Cannot put data from the offset bivector");
1835
1836 x_shift = cpl_vector_get(x_pos_offset,tel);
1837 y_shift = cpl_vector_get(y_pos_offset,tel);
1838
1839 /* In UV [m] */
1840 u_shift = (cfangle * x_shift - sfangle * y_shift)
1841 / scale;
1842 v_shift = (sfangle * x_shift + cfangle * y_shift)
1843 / scale;
1844 opd_pupil = -(u_shift * sobj_x + v_shift * sobj_y)
1846
1847 CPLCHECK_MSG("Cannot prepare data to be put in the ACQ PUPIL table");
1848
1849
1850 if (CHECK_BIT(is_it_bad,tel)||CHECK_BIT(is_it_bad,tel+4))
1851 {
1852 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1853
1854 /* increasing counters */
1855 if (CHECK_BIT(is_it_bad,tel )) n_bad_snr += 1;
1856 if (CHECK_BIT(is_it_bad,tel+4)) n_bad_distance += 1;
1857 u_shift = 0.0;
1858 v_shift = 0.0;
1859 opd_pupil = 0.0;
1860 cpl_msg_info(cpl_func,"Pupil image number%4lli is un-usable for tel %i",row,tel);
1861 }
1862 else
1863 {
1864 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel,
1865 16);
1866
1867 /* averaging positions */
1868 x_shift_sum += x_shift;
1869 y_shift_sum += y_shift;
1870 u_shift_sum += u_shift;
1871 v_shift_sum += v_shift;
1872
1873 /* increasing counter */
1874 n_sum += 1;
1875 }
1876
1877 cpl_table_set(acqcam_table, "PUPIL_X", row * GRAVI_SPOT_NTEL + tel,
1878 x_shift);
1879 cpl_table_set(acqcam_table, "PUPIL_Y", row * GRAVI_SPOT_NTEL + tel,
1880 y_shift);
1881 cpl_table_set(acqcam_table, "PUPIL_Z", row * GRAVI_SPOT_NTEL + tel,
1882 z_shift);
1883 cpl_table_set(acqcam_table, "PUPIL_R", row * GRAVI_SPOT_NTEL + tel,
1884 r_shift);
1885 cpl_table_set(acqcam_table, "PUPIL_U", row * GRAVI_SPOT_NTEL + tel,
1886 u_shift);
1887 cpl_table_set(acqcam_table, "PUPIL_V", row * GRAVI_SPOT_NTEL + tel,
1888 v_shift);
1889 cpl_table_set(acqcam_table, "PUPIL_W", row * GRAVI_SPOT_NTEL + tel,
1890 w_shift);
1891 cpl_table_set(acqcam_table, "OPD_PUPIL", row * GRAVI_SPOT_NTEL + tel,
1892 opd_pupil);
1893 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1894
1895 n_on++;
1896 }
1897 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1898 }
1899
1900 /* Add QC parameters */
1901 char qc_name[100];
1902
1903 sprintf(qc_name, "ESO QC ACQ FIELD%i NORTH_ANGLE", tel + 1);
1904 cpl_msg_info(cpl_func, "%s = %f", qc_name, northangle);
1905 cpl_propertylist_update_double(o_header, qc_name, northangle);
1906 cpl_propertylist_set_comment(o_header, qc_name,
1907 "[deg] y->x, predicted North direction on ACQ");
1908
1909 sprintf(qc_name, "ESO QC ACQ PUP%i FRAMES", tel + 1);
1910 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_sum);
1911 cpl_propertylist_update_double(o_header, qc_name, n_sum);
1912 cpl_propertylist_set_comment(o_header, qc_name,
1913 "Good ACQ pupil frames");
1914
1915 sprintf(qc_name, "ESO QC ACQ PUP%i BADSNR", tel + 1);
1916 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_bad_snr);
1917 cpl_propertylist_update_double(o_header, qc_name, n_bad_snr);
1918 cpl_propertylist_set_comment(o_header, qc_name,
1919 "Frames with low SNR");
1920
1921 sprintf(qc_name, "ESO QC ACQ PUP%i JUMP", tel + 1);
1922 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_bad_distance);
1923 cpl_propertylist_update_double(o_header, qc_name, n_bad_distance);
1924 cpl_propertylist_set_comment(o_header, qc_name,
1925 "Frames with position jump");
1926
1927 sprintf(qc_name, "ESO QC ACQ PUP%i SCALE", tel + 1);
1928 cpl_msg_info(cpl_func, "%s = %f", qc_name, scale);
1929 cpl_propertylist_update_double(o_header, qc_name, scale);
1930 cpl_propertylist_set_comment(o_header, qc_name,
1931 "[pix/m] diode scale on ACQ");
1932
1933 sprintf(qc_name, "ESO QC ACQ PUP%i FOCUS", tel + 1);
1934 cpl_msg_info(cpl_func, "%s = %f", qc_name, z_shift);
1935 cpl_propertylist_update_int(o_header, qc_name, z_shift);
1936 cpl_propertylist_set_comment(o_header, qc_name,
1937 "[pix] defocus of pupil plane");
1938
1939 /* to avoid division by zero */
1940 if (n_sum ==0) n_sum+=1;
1941
1942 sprintf(qc_name, "ESO QC ACQ PUP%i XPOS", tel + 1);
1943 cpl_msg_info(cpl_func, "%s = %f", qc_name, x_shift_sum/n_sum);
1944 cpl_propertylist_update_double(o_header, qc_name, x_shift_sum/n_sum);
1945 cpl_propertylist_set_comment(o_header, qc_name,
1946 "[pix] pupil x-shift in ACQ");
1947
1948 sprintf(qc_name, "ESO QC ACQ PUP%i YPOS", tel + 1);
1949 cpl_msg_info(cpl_func, "%s = %f", qc_name, y_shift_sum/n_sum);
1950 cpl_propertylist_update_double(o_header, qc_name, y_shift_sum/n_sum);
1951 cpl_propertylist_set_comment(o_header, qc_name,
1952 "[pix] pupil y-shift in ACQ");
1953
1954 sprintf(qc_name, "ESO QC ACQ PUP%i UPOS", tel + 1);
1955 cpl_msg_info(cpl_func, "%s = %f", qc_name, u_shift_sum/n_sum);
1956 cpl_propertylist_update_double(o_header, qc_name, u_shift_sum/n_sum);
1957 cpl_propertylist_set_comment(o_header, qc_name,
1958 "[m] pupil u-shift in ACQ");
1959
1960 sprintf(qc_name, "ESO QC ACQ PUP%i VPOS", tel + 1);
1961 cpl_msg_info(cpl_func, "%s = %f", qc_name, v_shift_sum/n_sum);
1962 cpl_propertylist_update_double(o_header, qc_name, v_shift_sum/n_sum);
1963 cpl_propertylist_set_comment(o_header, qc_name,
1964 "[m] pupil v-shift in ACQ");
1965
1966
1967 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1968
1969
1970
1971}
1972
1974return CPL_ERROR_NONE;
1975}
1976
1977
1978
1979/*----------------------------------------------------------------------------*/
1996/*----------------------------------------------------------------------------*/
1997
1998
1999
2000cpl_image * gravi_image_extract(cpl_image * image_in, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
2001 {
2002
2003 cpl_size nx = cpl_image_get_size_x (image_in);
2004 cpl_size ny = cpl_image_get_size_y (image_in);
2005 cpl_size llx_new,lly_new,urx_new,ury_new;
2006 cpl_size xpos, ypos;
2007
2008 cpl_size nx_new=urx+1-llx;
2009 cpl_size ny_new=ury+1-lly;
2010 cpl_image * image_out= cpl_image_new (nx_new, ny_new, cpl_image_get_type(image_in));
2011 cpl_image_fill_window (image_out, 1, 1, nx_new, ny_new, 0.0);
2012
2013
2014 if ((llx >= nx)||(lly >= ny)||(urx <= 0)||(ury <= 0))
2015 {
2016 cpl_msg_warning (cpl_func, "Cutting at x=(%lli,%lli) y=(%lli,%lli) is outside the window bondaries", llx, urx, lly, ury);
2017 }
2018 else
2019 {
2020 if (llx < 1)
2021 {
2022 llx_new=1;
2023 xpos=1-llx;
2024 } else
2025 {
2026 llx_new=llx;
2027 xpos=1;
2028 }
2029 if (lly < 1)
2030 {
2031 lly_new=1;
2032 ypos=1-lly;
2033 } else
2034 {
2035 lly_new=lly;
2036 ypos=1;
2037 }
2038
2039 if (urx > nx)
2040 urx_new=nx;
2041 else
2042 urx_new=urx;
2043 if (ury > ny)
2044 ury_new=ny;
2045 else
2046 ury_new=ury;
2047
2048 cpl_image * image_in_cut = cpl_image_extract(image_in,llx_new,lly_new,urx_new,ury_new);
2049
2050 cpl_image_copy (image_out, image_in_cut, xpos, ypos);
2051
2052 cpl_image_delete(image_in_cut);
2053 }
2054
2055 return image_out;
2056 }
2057
2058/*----------------------------------------------------------------------------*/
2063/*----------------------------------------------------------------------------*/
2064
2065
2067{
2068 double defocus;
2069 defocus = 0.4*(2.0*focus/(GRAVI_SPOT_NFOCUS-1) - 1.0)*(2.0*focus/(GRAVI_SPOT_NFOCUS-1) - 1.0)*(2.0*focus/(GRAVI_SPOT_NFOCUS-1) - 1.0);
2070
2071 return defocus;
2072}
2073
2074/*----------------------------------------------------------------------------*/
2091/*----------------------------------------------------------------------------*/
2092
2093cpl_error_code gravi_acqcam_field (cpl_image * mean_img,
2094 cpl_imagelist * acqcam_imglist,
2095 cpl_propertylist * header,
2096 cpl_table * acqcam_table,
2097 cpl_propertylist * o_header,
2098 gravi_data *static_param_data)
2099{
2101 cpl_ensure_code (mean_img, CPL_ERROR_NULL_INPUT);
2102 cpl_ensure_code (acqcam_imglist, CPL_ERROR_NULL_INPUT);
2103 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
2104 cpl_ensure_code (acqcam_table, CPL_ERROR_NULL_INPUT);
2105 cpl_ensure_code (o_header, CPL_ERROR_NULL_INPUT);
2106
2107 char qc_name[100];
2108 int ntel = 4;
2109
2110 int sts_mode = 0;
2111
2112 /* check the feed mode */
2113 if ( !strcmp(gravi_pfits_get_feed (header), "DUAL_STS") ) {
2114 sts_mode = DUAL_STS;
2115 }
2116 else{
2117 sts_mode = SINGLE_STS;
2118 }
2119
2120 /* Number of row */
2121 cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
2122
2123 /* MODE_ONAXIS or MODE_OFFAIS. We query once
2124 at the begining of the function */
2125 int axis_mode = gravi_pfits_get_axis (header);
2126
2127 /* Create columns */
2128 gravi_table_new_column (acqcam_table, "FIELD_SC_X", "pix", CPL_TYPE_DOUBLE);
2129 gravi_table_new_column (acqcam_table, "FIELD_SC_Y", "pix", CPL_TYPE_DOUBLE);
2130 gravi_table_new_column (acqcam_table, "FIELD_SC_XERR", "pix", CPL_TYPE_DOUBLE);
2131 gravi_table_new_column (acqcam_table, "FIELD_SC_YERR", "pix", CPL_TYPE_DOUBLE);
2132 gravi_table_new_column (acqcam_table, "FIELD_FT_X", "pix", CPL_TYPE_DOUBLE);
2133 gravi_table_new_column (acqcam_table, "FIELD_FT_Y", "pix", CPL_TYPE_DOUBLE);
2134 gravi_table_new_column (acqcam_table, "FIELD_FT_XERR", "pix", CPL_TYPE_DOUBLE);
2135 gravi_table_new_column (acqcam_table, "FIELD_FT_YERR", "pix", CPL_TYPE_DOUBLE);
2136 gravi_table_new_column (acqcam_table, "FIELD_SCALE", "mas/pix", CPL_TYPE_DOUBLE);
2137 gravi_table_new_column (acqcam_table, "FIELD_SCALEERR", "mas/pix", CPL_TYPE_DOUBLE);
2138 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DX", "pix", CPL_TYPE_DOUBLE);
2139 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DY", "pix", CPL_TYPE_DOUBLE);
2140 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DXERR", "pix", CPL_TYPE_DOUBLE);
2141 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DYERR", "pix", CPL_TYPE_DOUBLE);
2142 gravi_table_new_column (acqcam_table, "FIELD_STREHL", "ratio", CPL_TYPE_DOUBLE);
2143
2144 /* ----------------- START EKW 26/11/2018 read constant parameter from calibration file */
2145 /* Position of roof center on full frame */
2146 //double roof_x[] = {274.4, 787.1, 1236.1, 1673.4};
2147 //double roof_y[] = {242.3, 247.7, 225.8, 235.6};
2148 double *roof_x = NULL;
2149 double *roof_y = NULL;
2150
2151 /* Position of single-field spot on full frame */
2152 //double spot_x[] = {289. , 798.2, 1245.5, 1696.};
2153 //double spot_y[] = {186.5, 187.5, 172.5, 178.};
2154 double *spot_x = NULL;
2155 double *spot_y = NULL;
2156
2157 /* Default position angle of roof */
2158 //double roof_pos[] = {38.49, 38.54, 38.76, 39.80};
2159 double *roof_pos = NULL;
2160
2161 cpl_table * roof_pos_table = gravi_data_get_table (static_param_data, "ROOFPOS");
2162 CPLCHECK_MSG ("STATIC_PARAM not available in the SOF. It is mandatory for acqcam reduction.");
2163
2164 if ( cpl_table_has_column(roof_pos_table , "roof_x") ) {
2165 roof_x= cpl_table_get_data_double (roof_pos_table, "roof_x");
2166 cpl_msg_info(cpl_func,"roof_x [0] : %e", roof_x[0] );
2167 cpl_msg_info(cpl_func,"roof_x [1] : %e", roof_x[1] );
2168 cpl_msg_info(cpl_func,"roof_x [2] : %e", roof_x[2] );
2169 cpl_msg_info(cpl_func,"roof_x [3] : %e", roof_x[3] );
2170 }
2171 else {
2172 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_x ");
2173 }
2174
2175 if ( cpl_table_has_column(roof_pos_table , "roof_y") ) {
2176 roof_y= cpl_table_get_data_double (roof_pos_table, "roof_y");
2177 cpl_msg_info(cpl_func,"roof_y [0] : %e", roof_y[0] );
2178 cpl_msg_info(cpl_func,"roof_y [1] : %e", roof_y[1] );
2179 cpl_msg_info(cpl_func,"roof_y [2] : %e", roof_y[2] );
2180 cpl_msg_info(cpl_func,"roof_y [3] : %e", roof_y[3] );
2181 }
2182 else {
2183 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_y ");
2184 }
2185
2186 if ( cpl_table_has_column(roof_pos_table , "spot_x") ) {
2187 spot_x= cpl_table_get_data_double (roof_pos_table, "spot_x");
2188 cpl_msg_info(cpl_func,"spot_x [0] : %e", spot_x[0] );
2189 cpl_msg_info(cpl_func,"spot_x [1] : %e", spot_x[1] );
2190 cpl_msg_info(cpl_func,"spot_x [2] : %e", spot_x[2] );
2191 cpl_msg_info(cpl_func,"spot_x [3] : %e", spot_x[3] );
2192 }
2193 else {
2194 cpl_msg_warning(cpl_func,"Cannot get the default values for spot_x ");
2195 }
2196
2197 if ( cpl_table_has_column(roof_pos_table , "spot_y") ) {
2198 spot_y= cpl_table_get_data_double (roof_pos_table, "spot_y");
2199 cpl_msg_info(cpl_func,"spot_y [0] : %e", spot_y[0] );
2200 cpl_msg_info(cpl_func,"spot_y [1] : %e", spot_y[1] );
2201 cpl_msg_info(cpl_func,"spot_y [2] : %e", spot_y[2] );
2202 cpl_msg_info(cpl_func,"spot_y [3] : %e", spot_y[3] );
2203 }
2204 else {
2205 cpl_msg_warning(cpl_func,"Cannot get the default values for spot_y ");
2206 }
2207
2208 if ( cpl_table_has_column(roof_pos_table , "roof_pos") ) {
2209 roof_pos= cpl_table_get_data_double (roof_pos_table, "roof_pos");
2210 cpl_msg_info(cpl_func,"roof_pos [0] : %e", roof_pos[0] );
2211 cpl_msg_info(cpl_func,"roof_pos [1] : %e", roof_pos[1] );
2212 cpl_msg_info(cpl_func,"roof_pos [2] : %e", roof_pos[2] );
2213 cpl_msg_info(cpl_func,"roof_pos [3] : %e", roof_pos[3] );
2214 }
2215 else {
2216 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_pos ");
2217 }
2218 /* ------------------ END EKW 26/11/2018 read constant parameter from calibration file */
2219
2220
2221
2222 /* Static acqcam field calibration error */
2223 /* calibration from Frank, March 2017, using Marcel
2224 double calib_dx[] = {-0.3483, -1.0251, -0.5432, -0.2024} ;
2225 double calib_dy[] = { 0.3089, -0.5221, -0.2686, -0.3843} ; */
2226 /* calibration from Julien, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-10
2227 double calib_dx[] = {-0.61, -1.69, -0.97, 0.00} ;
2228 double calib_dy[] = { 0.53, -1.07, -0.49, -0.47} ; */
2229 /* calibration from Oli, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-19/
2230 double calib_dx[] = {-0.63, -1.67, -0.97, 0.02} ;
2231 double calib_dy[] = { 0.55, -1.07, -0.49, -0.45} ; */
2232 /* FE: removing obsolete calib_dx, because now handled by ICS, or by
2233 correcting fitsheader retrospectively */
2234 double calib_dx[] = {0,0,0,0} ;
2235 double calib_dy[] = {0,0,0,0} ;
2236
2237 /* If sub-windowing, we read the sub-window size */
2238 cpl_size nsx = 512;
2239 /* EKW 10/01/2019 cpl_size nsy = 512; */
2240 if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
2241 nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
2242 /* EKW 10/01/2019 nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY"); */
2243 }
2244
2245 /* Compute separation */
2246 double sobj_x = gravi_pfits_get_sobj_x (header);
2247 double sobj_y = gravi_pfits_get_sobj_y (header);
2248
2249 /* Accumulated mapping offsets */
2250 double sobj_offx=0., sobj_offy=0, sobj_drho=0., sobj_dth=0.;
2251 if (cpl_propertylist_has(header, "ESO INS SOBJ OFFX")) {
2252 /* mapping / mosaicing blind offsets from otiginal acquisition
2253 * to current position */
2254 sobj_offx = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFX");
2255 sobj_offy = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFY");
2256
2257 /* distance from acquired SC position to current SC position in mas */
2258 sobj_drho = sqrt(sobj_offx*sobj_offx+sobj_offy*sobj_offy);
2259
2260 /* position angle on sky from acquired SC position to */
2261 /* current SC position, neglecting anamorphism variations */
2262 sobj_dth = atan2(sobj_offx, sobj_offy)*180./M_PI;
2263 }
2264 CPLCHECK_MSG ("Cannot get separation");
2265
2266 /* Recover position of originally acquired star, before any blind offset */
2267 double dx_in = sobj_x - sobj_offx;
2268 double dy_in = sobj_y - sobj_offy;
2269 double rho_in = sqrt(dx_in*dx_in + dy_in*dy_in);
2270
2271 /* Force zero separation for SINGLE mode */
2272 if ( gravi_pfits_get_mode (header) == MODE_SINGLE) rho_in = 0.;
2273 CPLCHECK_MSG ("Error reading header information");
2274
2275 /* Loop on tel */
2276 for (int tel = 0; tel < ntel; tel++) {
2277 char name[90];
2278 double rp=roof_pos[tel]; // default value of roof position angle
2279 cpl_size sx=tel*512+1, sy=1;
2280 cpl_msg_info (cpl_func, "Compute field position for beam %i", tel+1);
2281
2282 /* Read roof position angle */
2283 sprintf (name, "ESO INS DROTOFF%d", tel + 1);
2284 if ( cpl_propertylist_has (header, name) ) {
2285 rp = cpl_propertylist_get_double(header, name);
2286 CPLCHECK ("Cannot get rotation");
2287 }
2288 else {
2289 cpl_msg_info (cpl_func, "%s not in header, use %f", name, rp);
2290 }
2291
2292 /* Approx. position angle of the binary, left from top */
2293 double approx_PA = 270.-rp;
2294
2295 /* Get the telescope name and ID */
2296 const char * telname = gravi_conf_get_telname (tel, header);
2297
2298 /* Check telescope name */
2299 if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
2300 cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
2301
2302 /* Hardcoded approx. plate-scale in mas/pix. */
2303 double scale = 0.0;
2304 char telid[128];
2305 if (telname[0] == 'U') {
2306 /* scale = 18.; */
2307 scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO PLATE SCALE UT");
2308 cpl_msg_info (cpl_func,"PLATE SCALE UT is : %e",scale);
2309 } else if (telname[0] == 'A') {
2310 sprintf(telid, "ESO PLATE SCALE AT%d", tel+1);
2311 scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), telid);
2312 cpl_msg_info (cpl_func,"PLATE SCALE AT%d is : %e",tel+1, scale);
2313 /*
2314 if (tel == 0) {
2315 scale = 76.8;
2316 } else if (tel == 1) {
2317 scale = 78.0;
2318 } else if (tel == 2) {
2319 scale = 77.0;
2320 } else if (tel == 3) {
2321 scale = 84.6;
2322 }
2323 */
2324 }
2325
2326 /* Position of the fibres */
2327 double fiber_xft=0.;
2328 double fiber_yft=0.;
2329 double fiber_xsc=0.;
2330 double fiber_ysc=0.;
2331 sprintf (name, "ESO ACQ FIBER FT%dX", tel + 1);
2332 if (cpl_propertylist_has(header, name)) {
2333 fiber_xft=cpl_propertylist_get_double(header, name);
2334 sprintf (name, "ESO ACQ FIBER FT%dY", tel + 1);
2335 fiber_yft=cpl_propertylist_get_double(header, name);
2336 sprintf (name, "ESO ACQ FIBER SC%dX", tel + 1);
2337 fiber_xsc=cpl_propertylist_get_double(header, name);
2338 sprintf (name, "ESO ACQ FIBER SC%dY", tel + 1);
2339 fiber_ysc=cpl_propertylist_get_double(header, name);
2340 }
2341 double fiber_ft_sc_x=fiber_xsc-fiber_xft;
2342 double fiber_ft_sc_y=fiber_ysc-fiber_yft;
2343
2344 /* Get the North position angle on the camera */
2345 double northangle = gravi_pfits_get_northangle_acqcam (header, tel);
2346 CPLCHECK ("Cannot get rotation");
2347
2348 /* Mapping/mosaicing offset on acq cam axes, in mas, */
2349 /* neglecting amnamorphism variations */
2350 double sobj_offx_cam = sobj_drho * sin((northangle+sobj_dth)/180.*M_PI);
2351 double sobj_offy_cam = sobj_drho * cos((northangle+sobj_dth)/180.*M_PI);
2352
2353 /* If sub-windowing, we read the sub-window start for field */
2354 if ( nsx != 512 ) {
2355 sprintf (name, "ESO DET1 FRAM%d STRX", tel + 1);
2356 sx = cpl_propertylist_get_int (header, name);
2357
2358 sprintf (name, "ESO DET1 FRAM%d STRY", tel + 1);
2359 sy = cpl_propertylist_get_int (header, name);
2360
2361 CPLCHECK_MSG ("Cannot get sub-windowing parameters");
2362 }
2363
2364 cpl_msg_debug (cpl_func,"sub-window field %i sx= %lld sy = %lld", tel, sx, sy);
2365
2366 /*--------------------------------------------------*/
2367 /* Guess position of SC and FT targes in mean image */
2368 /*--------------------------------------------------*/
2369
2370 double xFT, yFT, xSC, ySC;
2371
2372 if (axis_mode == MODE_ONAXIS) {
2373 /* TODO: close dual-field */
2374 /* Single-field case */
2375 /* Simply shift the best spot from full frame to cut-out */
2376 xFT = spot_x[tel] - sx + nsx*tel + 1;
2377 yFT = spot_y[tel] - sy + 1;
2378 xSC = xFT;
2379 ySC = yFT;
2380 } else if (sts_mode == DUAL_STS) {
2381 xFT = fiber_xft - sx + nsx*tel + 1;
2382 yFT = fiber_yft - sy + 1;
2383 xSC = fiber_xsc - sx + nsx*tel + 1;
2384 ySC = fiber_ysc - sy + 1;
2385 } else {
2386 /* Pixel position of roof center on cut-out frame */
2387 /* Shift from full frame to cut-out */
2388 double cutout_roof_x = roof_x[tel] - sx + nsx*tel + 1;
2389 double cutout_roof_y = roof_y[tel] - sy + 1;
2390
2391 cpl_msg_info (cpl_func, " ROOF_X = %f", cutout_roof_x);
2392 cpl_msg_info (cpl_func, " ROOF_Y = %f", cutout_roof_y);
2393
2394 /* Approx pixel offset from SC to FT, divided by 2 */
2395 double approx_dx=0.5*rho_in*sin(approx_PA*M_PI/180.)/scale;
2396 double approx_dy=0.5*rho_in*cos(approx_PA*M_PI/180.)/scale;
2397
2398 /* Expected position of the two stars */
2399 xFT = cutout_roof_x - approx_dx ;
2400 yFT = cutout_roof_y - approx_dy ;
2401 xSC = cutout_roof_x + approx_dx ;
2402 ySC = cutout_roof_y + approx_dy ;
2403 }
2404
2405 cpl_msg_info (cpl_func, "guess SC_X = %f", xSC);
2406 cpl_msg_info (cpl_func, "guess SC_Y = %f", ySC);
2407 cpl_msg_info (cpl_func, "guess FT_X = %f", xFT);
2408 cpl_msg_info (cpl_func, "guess FT_Y = %f", yFT);
2409
2410 /*------------------------------------------*/
2411 /* Measure SC target position in mean image */
2412 /*------------------------------------------*/
2413
2414 /* Box size */
2415 /* Optimal size has been determined empirically. */
2416 /* Too small value (15) will sometimes miss even a bright target. */
2417 /* Too large value (25) will often pick a wrong star (or backround noise) */
2418 cpl_size size = 20;
2419
2420 double xSCguess=xSC, ySCguess=ySC, xFTguess=xFT, yFTguess=yFT;
2421 double ex, ey;
2422 double qc_val=0.;
2423
2424 gravi_acq_fit_gaussian (mean_img, &xSC, &ySC, &ex, &ey, size);
2425 CPLCHECK_MSG("Error fitting SC");
2426
2427 /*--------------------------------------------------------*/
2428 /* Add QC parameters for SC target position in mean image */
2429 /*--------------------------------------------------------*/
2430
2431 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X", tel+1);
2432 if (xSC==0.) {
2433 /* Fitting failed: put QC to 0., reset xSC to gues value */
2434 qc_val = 0.;
2435 xSC = xSCguess;
2436 } else {
2437 /* Fiting succeeded: shift into full frame */
2438 qc_val = xSC;
2439 }
2440 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2441 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2442 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2443
2444 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y", tel+1);
2445 if (ySC==0.) {
2446 /* Fitting failed: put QC to 0., reset ySC to gues value */
2447 qc_val = 0.;
2448 ySC = ySCguess;
2449 } else {
2450 /* Fiting succeeded: shift into full frame */
2451 qc_val = ySC;
2452 }
2453 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2454 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2455 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2456
2457 /*------------------------------------------*/
2458 /* Measure FT target position in mean image */
2459 /*------------------------------------------*/
2460
2461 if (axis_mode != MODE_ONAXIS) {
2462 /* Detec FT */
2463 gravi_acq_fit_gaussian (mean_img, &xFT, &yFT, &ex, &ey, size);
2464 CPLCHECK_MSG("Error fitting FT");
2465 } else {
2466 xFT=xSC;
2467 yFT=ySC;
2468 }
2469
2470 /*--------------------------------------------------------*/
2471 /* Add QC parameters for FT target position in mean image */
2472 /*--------------------------------------------------------*/
2473
2474 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X", tel+1);
2475 if (xFT==0.) {
2476 /* Fitting failed: put QC to 0., reset xFT to gues value */
2477 qc_val = 0.;
2478 xFT = xFTguess;
2479 } else {
2480 /* Fiting succeeded: shift into full frame */
2481 qc_val = xFT;
2482 }
2483 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2484 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2485 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2486
2487 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y", tel+1);
2488 if (yFT==0.) {
2489 /* Fitting failed: put QC to 0., reset yFT to gues value */
2490 qc_val = 0.;
2491 yFT = yFTguess;
2492 } else {
2493 /* Fiting succeeded: shift into full frame */
2494 qc_val = yFT;
2495 }
2496 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2497 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2498 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2499
2500 /*---------------------------------------*/
2501 /* If in dual field, measure plate scale */
2502 /*---------------------------------------*/
2503
2504 if (axis_mode != MODE_ONAXIS) {
2505 sprintf (qc_name, "ESO QC ACQ FIELD%i SCALE", tel+1);
2506 double sep = sqrt((ySC-yFT)*(ySC-yFT)+(xSC-xFT)*(xSC-xFT));
2507 /* FE 2019-10-31 in case separation could not be measured, i.e. 0, then use default scale */
2508 double pscale = scale;
2509 if (sts_mode != DUAL_STS) pscale = sep ? rho_in/sep : scale;
2510 /* double pscale = sep ? rho_in/sep : 0.; */
2511 qc_val=pscale;
2512 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2513 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2514 cpl_propertylist_set_comment (o_header, qc_name,
2515 "[mas/pixel] plate-scale in the "
2516 "FT-SC direction");
2517 }
2518
2519 /*------------------------------------------------*/
2520 /* If in dual field, measure fibre position error */
2521 /*------------------------------------------------*/
2522
2523 if (axis_mode != MODE_ONAXIS) {
2524 /* Error in SC fibre positioning */
2525 /* The three terms are */
2526 /* - offset from FT target as detected to original SC */
2527 /* target as detected; */
2528 /* - blind offset command from mapping template projected */
2529 /* on acqisition camera; */
2530 /* - offset from FT fiber to SC fiber. */
2531 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DX", tel+1);
2532 qc_val = 0;
2533 if (scale) {
2534 qc_val = (xSC-xFT) + sobj_offx_cam/scale - fiber_ft_sc_x - calib_dx[tel];
2535 }
2536 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2537 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2538 cpl_propertylist_set_comment (o_header, qc_name,
2539 "[pixel] dx from SC fiber to "
2540 "SC object");
2541
2542 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DY", tel+1);
2543 qc_val = 0;
2544 if (scale) {
2545 qc_val = (ySC-yFT) + sobj_offy_cam/scale - fiber_ft_sc_y - calib_dy[tel];
2546 }
2547 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2548 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2549 cpl_propertylist_set_comment (o_header, qc_name,
2550 "[pixel] dx from SC fiber to "
2551 "SC object");
2552 }
2553
2554 /*-----------------------------*/
2555 /* Average Strehl computation */
2556 /*-----------------------------*/
2557
2558 double max_on_average, strehl_on_average;
2559
2560 /* Measure strehl and maximum on averaged image */
2561 gravi_acq_measure_strehl(mean_img, xFT, yFT, scale, &strehl_on_average, header);
2562
2563 /*
2564 * JIRA PIPE-9123 : Pipeline faile with Strehl=NaN
2565 * For 2017 AT data the spots jump up by 30 pixel and therefore are not
2566 * catched in the subwindow around roof_pos , cut_out +/- 50 pixel
2567 */
2568 if (cpl_error_get_code() != CPL_ERROR_NONE ) {
2569 cpl_msg_info (cpl_func, "WARNING Filling STREHL due to NaN");
2570 cpl_error_reset();
2571 strehl_on_average = 0.0;
2572 }
2573
2574 gravi_acq_measure_max(mean_img, xFT, yFT, 15, &max_on_average);
2575
2576 /* Update Strehl QC */
2577 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL", tel+1);
2578 cpl_msg_info (cpl_func, "%s = %f", qc_name, strehl_on_average);
2579 cpl_propertylist_update_double (o_header, qc_name, strehl_on_average);
2580 cpl_propertylist_set_comment (o_header, qc_name, "Average Strehl from stacked AcqCam images");
2581
2582
2583 /* Adding Strehl for SC channel in wide mode (PIPE-9913) */
2584 if (sts_mode == DUAL_STS) {
2585 gravi_acq_measure_strehl(mean_img, xSC, ySC, scale, &strehl_on_average, header);
2586
2587 /*
2588 * JIRA PIPE-9123 : Pipeline faile with Strehl=NaN
2589 * For 2017 AT data the spots jump up by 30 pixel and therefore are not
2590 * catched in the subwindow around roof_pos , cut_out +/- 50 pixel
2591 */
2592 if (cpl_error_get_code() != CPL_ERROR_NONE ) {
2593 cpl_msg_info (cpl_func, "WARNING Filling STREHL SC due to NaN");
2594 cpl_error_reset();
2595 strehl_on_average = 0.0;
2596 }
2597
2598 /* Update Strehl QC */
2599 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHLSC", tel+1);
2600 cpl_msg_info (cpl_func, "%s = %f", qc_name, strehl_on_average);
2601 cpl_propertylist_update_double (o_header, qc_name, strehl_on_average);
2602 cpl_propertylist_set_comment (o_header, qc_name, "Average Strehl from stacked AcqCam images");
2603 }
2604 /*----------------------------*/
2605 /* Now process frame by frame */
2606 /*----------------------------*/
2607
2608 /* Loop on all images */
2609 for (cpl_size row = 0; row < nrow; row++) {
2610 if (row %10 == 0 || row == (nrow-1))
2611 cpl_msg_info_overwritable (cpl_func, "Fit image %lld over %lld", row+1, nrow);
2612
2613 /* Get data */
2614 cpl_image * img = cpl_imagelist_get (acqcam_imglist, row);
2615 CPLCHECK_MSG("Error getting image");
2616
2617 /*-------------------------------------------------*/
2618 /* SC target position computation of current frame */
2619 /*-------------------------------------------------*/
2620 double xsc = xSC, ysc = ySC, exsc=0., eysc=0.;
2621 gravi_acq_fit_gaussian (img, &xsc, &ysc, &exsc, &eysc, size);
2622 CPLCHECK_MSG("Error fitting SC");
2623 /* Shift back positions to full frame */
2624 if (xsc != 0.) xsc += sx - 1 - nsx*tel;
2625 if (ysc != 0.) ysc += sy - 1;
2626 cpl_table_set (acqcam_table, "FIELD_SC_X", row*ntel+tel, xsc);
2627 cpl_table_set (acqcam_table, "FIELD_SC_Y", row*ntel+tel, ysc);
2628 cpl_table_set (acqcam_table, "FIELD_SC_XERR", row*ntel+tel, exsc);
2629 cpl_table_set (acqcam_table, "FIELD_SC_YERR", row*ntel+tel, eysc);
2630 CPLCHECK_MSG("Error setting SC columns");
2631
2632 /*-------------------------------------------------*/
2633 /* FT target position computation of current frame */
2634 /*-------------------------------------------------*/
2635 double xft = xFT, yft = yFT, exft=0., eyft=0.;
2636 if (axis_mode != MODE_ONAXIS) {
2637 gravi_acq_fit_gaussian (img, &xft, &yft, &exft, &eyft, size);
2638 CPLCHECK_MSG("Error fitting FT");
2639 /* Shift back positions to full frame */
2640 if (xft != 0.) xft += sx - 1 - nsx*tel;
2641 if (yft != 0.) yft += sy - 1;
2642 } else {
2643 xft=xsc;
2644 yft=ysc;
2645 exft=exsc;
2646 eyft=eysc;
2647 }
2648 cpl_table_set (acqcam_table, "FIELD_FT_X", row*ntel+tel, xft);
2649 cpl_table_set (acqcam_table, "FIELD_FT_Y", row*ntel+tel, yft);
2650 cpl_table_set (acqcam_table, "FIELD_FT_XERR", row*ntel+tel, exft);
2651 cpl_table_set (acqcam_table, "FIELD_FT_YERR", row*ntel+tel, eyft);
2652 CPLCHECK_MSG("Error setting FT column");
2653
2654 /*------------------------------------------*/
2655 /* Plate scale computation of current frame */
2656 /*------------------------------------------*/
2657 double ft_sc_x = xsc - xft;
2658 double ft_sc_y = ysc - yft;
2659 double eft_sc_x = sqrt(exsc*exsc+exft*exft);
2660 double eft_sc_y = sqrt(eysc*eysc+eyft*eyft);
2661 double sep = sqrt(ft_sc_x*ft_sc_x+ft_sc_y*ft_sc_y);
2662 double pscale = scale;
2663 if (sts_mode != DUAL_STS) pscale = sep ? rho_in/sep : 0.;
2664 double escale = 0.;
2665 if (sep) {
2666 escale = rho_in/(sep*sep*sep)*
2667 (ft_sc_x*eft_sc_x+ft_sc_y*eft_sc_y) ;
2668 }
2669 cpl_table_set (acqcam_table, "FIELD_SCALE", row*ntel+tel, pscale);
2670 cpl_table_set (acqcam_table, "FIELD_SCALEERR", row*ntel+tel, escale);
2671
2672 /*---------------------------------------------*/
2673 /* Error in fibre positioning of current frame */
2674 /*---------------------------------------------*/
2675 /* The three terms are */
2676 /* - offset from FT target as detected to original SC */
2677 /* target as detected; */
2678 /* - blind offset command from mapping template projected */
2679 /* on acqisition camera; */
2680 /* - offset from FT fiber to SC fiber. */
2681 double corrx=0, corry=0., ecorrx=0., ecorry=0.;
2682 if (pscale) {
2683 corrx = ft_sc_x + sobj_offx_cam/pscale - fiber_ft_sc_x - calib_dx[tel];
2684 corry = ft_sc_y + sobj_offy_cam/pscale - fiber_ft_sc_y - calib_dy[tel];
2685 double tmp = escale/(pscale*pscale);
2686 tmp *= tmp;
2687 ecorrx = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offx_cam*sobj_offx_cam*tmp);
2688 ecorry = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offy_cam*sobj_offy_cam*tmp);
2689 }
2690 cpl_table_set (acqcam_table, "FIELD_FIBER_DX", row*ntel+tel, corrx);
2691 cpl_table_set (acqcam_table, "FIELD_FIBER_DY", row*ntel+tel, corry);
2692 cpl_table_set (acqcam_table, "FIELD_FIBER_DXERR", row*ntel+tel, ecorrx);
2693 cpl_table_set (acqcam_table, "FIELD_FIBER_DYERR", row*ntel+tel, ecorry);
2694
2695 /*-------------------------------------*/
2696 /* Strehl computation of current frame */
2697 /*-------------------------------------*/
2698 double max_on_frame;
2699 gravi_acq_measure_max(img, xFT, yFT, 15, &max_on_frame);
2700 cpl_table_set (acqcam_table, "FIELD_STREHL", row*ntel+tel, strehl_on_average*(max_on_frame/max_on_average) );
2701
2702 } /* End loop on images */
2703
2704 /* Add some QC */
2705 double sc_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_SC_X", tel, ntel);
2706 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X STD", tel+1);
2707 cpl_propertylist_update_double (o_header, qc_name, sc_std_x);
2708 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
2709
2710 double sc_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_SC_Y", tel, ntel);
2711 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y STD", tel+1);
2712 cpl_propertylist_update_double (o_header, qc_name, sc_std_y);
2713 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
2714
2715 double ft_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_FT_X", tel, ntel);
2716 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X STD", tel+1);
2717 cpl_propertylist_update_double (o_header, qc_name, ft_std_x);
2718 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
2719
2720 double ft_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_FT_Y", tel, ntel);
2721 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y STD", tel+1);
2722 cpl_propertylist_update_double (o_header, qc_name, ft_std_y);
2723 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
2724
2725 double strehl_std = gravi_table_get_column_std (acqcam_table, "FIELD_STREHL", tel, ntel);
2726 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL STD", tel+1);
2727 cpl_propertylist_update_double (o_header, qc_name, strehl_std);
2728 cpl_propertylist_set_comment (o_header, qc_name, "Std of FT strehl");
2729
2730 } /* End loop on tel */
2731
2733 return CPL_ERROR_NONE;
2734}
2735
2736/*----------------------------------------------------------------------------*/
2760/*----------------------------------------------------------------------------*/
2761
2762cpl_error_code gravi_reduce_acqcam (gravi_data * output_data,
2763 gravi_data * input_data,
2764 gravi_data * sky_data,
2765 gravi_data * dark_data,
2766 gravi_data * static_param_data,
2767 int saveAcqTable)
2768{
2770 cpl_ensure_code (output_data, CPL_ERROR_NULL_INPUT);
2771 cpl_ensure_code (input_data, CPL_ERROR_NULL_INPUT);
2772
2773 /* Check if extension exist */
2775 gravi_msg_warning (cpl_func, "Cannot reduce the ACQCAM, not data");
2776 return CPL_ERROR_NONE;
2777 }
2778
2779 /* Get the data and header */
2780 cpl_propertylist * header, * o_header;
2781 header = gravi_data_get_header (input_data);
2782 o_header = gravi_data_get_header (output_data);
2783
2784 /* Check if ISS data present (to avoid crash) */
2785 if (!gravi_conf_get_telname (0, header)) {
2786 gravi_msg_warning (cpl_func, "Cannot reduce the ACQCAM, no ISS keywords");
2787 return CPL_ERROR_NONE;
2788 }
2789
2790 cpl_imagelist * acqcam_imglist;
2791 cpl_imagelist * acqcam_imglist_v2;
2792 acqcam_imglist = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2793 acqcam_imglist_v2 = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2794 CPLCHECK_MSG ("Cannot get data or header");
2795
2796 /* Build the table */
2797 cpl_size ntel = 4;
2798 cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
2799
2800 cpl_table * acqcam_table;
2801 acqcam_table = cpl_table_new (nrow * ntel);
2802
2803 /* Time column shall contain the time from PCR.ACQ.START in [us] */
2804 cpl_table_new_column (acqcam_table, "TIME", CPL_TYPE_INT);
2805 cpl_table_set_column_unit (acqcam_table, "TIME", "us");
2806
2807 /* Loop on DIT in cube to fill the TIME column
2808 * same value for all 4 beams*/
2809 for (cpl_size row = 0; row < nrow; row++) {
2810 double time = gravi_pfits_get_time_acqcam (header, row);
2811 for (int tel = 0; tel < ntel; tel ++)
2812 cpl_table_set (acqcam_table, "TIME", row*ntel+tel, time);
2813 }
2814
2815 /* Compute mean image */
2816 cpl_image * mean_img = cpl_imagelist_collapse_create (acqcam_imglist);
2817 cpl_image * mean_subtracted_img = NULL;
2818
2819 if (sky_data) {
2820 mean_subtracted_img = gravi_data_get_img(sky_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2821 cpl_image_subtract(mean_subtracted_img, mean_img);
2822 cpl_image_multiply_scalar(mean_subtracted_img, -1.0);
2823 } else if (dark_data) {
2824 mean_subtracted_img = gravi_data_get_img(dark_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2825 cpl_msg_warning(cpl_func, "No SKY, using DARK for acquisition field subtraction");
2826 cpl_image_subtract(mean_subtracted_img, mean_img);
2827 cpl_image_multiply_scalar(mean_subtracted_img, -1.0);
2828 } else {
2829 cpl_msg_warning(cpl_func, "No SKY or DARK, skipping acquisition field subtraction");
2830 mean_subtracted_img = mean_img;
2831 }
2832
2833 /* Compute FIELD columns */
2834 gravi_acqcam_field (mean_subtracted_img, acqcam_imglist, header,
2835 acqcam_table, o_header, static_param_data);
2836 CPLCHECK_MSG ("Cannot reduce acquisition field images");
2837
2838 /* Compute the PUPIL and save output only if requested */
2839 if (saveAcqTable)
2840 {
2841 /* Compute PUPIL columns with algorithm 2.0 */
2842 gravi_acqcam_pupil_v2 (mean_img, acqcam_imglist_v2, header,
2843 acqcam_table, o_header, static_param_data);
2844 CPLCHECK_MSG ("Cannot reduce pupil images");
2845
2846 /* Add this output table in the gravi_data */
2847 cpl_propertylist * plist_acq_cam = cpl_propertylist_new ();
2848 cpl_propertylist_update_string (plist_acq_cam, "INSNAME", INSNAME_ACQ);
2849 gravi_data_add_img (output_data, plist_acq_cam, GRAVI_IMAGING_DATA_ACQ_EXT, mean_img);
2850 CPLCHECK_MSG ("Cannot add acqcam_table in data");
2851
2852 gravi_data_add_table (output_data, NULL, GRAVI_OI_VIS_ACQ_EXT, acqcam_table);
2853 CPLCHECK_MSG ("Cannot add acqcam_table in data");
2854 }
2855 else
2856 {
2857 FREE (cpl_table_delete, acqcam_table);
2858 FREE (cpl_image_delete, mean_img);
2859 }
2860
2862 return CPL_ERROR_NONE;
2863}
2864
2865/*----------------------------------------------------------------------------*/
2870/*----------------------------------------------------------------------------*/
2871
2872double gravi_acqcam_z2meter (double PositionPixels, gravi_data * static_param_data)
2873{
2874 double f_PT ; /* = 14e-3; */ /* pupil tracker lenslet FL*/
2875 double f_lens ; /* = 467e-3; */ /* folding optics lens FL */
2876 double Llambda ; /* = 1.2e-6; */ /* laser diode wavelength */
2877 double D_beam ; /* = 18e-3; */ /* meter */
2878 double D_pixel ; /* = 18e-6; */
2879 double D_AT ; /* = 1.8; */ /* m */
2880 double D_lenslet ; /* = 2 * 1.015e-3; */
2881
2882 f_PT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS fPT");
2883 cpl_msg_debug (cpl_func,"ESO LENS fPT is : %e", f_PT);
2884
2885 f_lens = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS flens");
2886 cpl_msg_debug (cpl_func,"ESO LENS f_lens is : %e", f_lens);
2887
2888 Llambda = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Llambda");
2889 cpl_msg_debug (cpl_func,"ESO LENS Llambda is : %e", Llambda);
2890
2891 D_beam = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dbeam");
2892 cpl_msg_debug (cpl_func,"ESO LENS D_beam is : %e", D_beam);
2893
2894 D_pixel = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dpixel");
2895 cpl_msg_debug (cpl_func,"ESO LENS D_pixel is : %e", D_pixel);
2896
2897 D_AT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS DAT");
2898 cpl_msg_debug (cpl_func,"ESO LENS D_AT is : %e", D_AT);
2899
2900 D_lenslet = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dlenslet");
2901 cpl_msg_debug (cpl_func,"ESO LENS D_lenslet is : %e", D_lenslet);
2902
2903
2904 double longDef;
2905 longDef = 8 * (f_PT / D_lenslet) * (f_PT / D_lenslet) * 3.5 * D_pixel *
2906 D_beam / (f_PT * D_lenslet) * Llambda / CPL_MATH_2PI * PositionPixels;
2907
2908 return f_lens * f_lens * longDef / (f_PT + longDef) / f_PT * (D_AT / D_lenslet);
2909}
2910
2911
2912/*----------------------------------------------------------------------------*/
2923/*----------------------------------------------------------------------------*/
2924
2925cpl_error_code gravi_image_fft_correlate (cpl_image *ia, cpl_image *ib, cpl_size *xd, cpl_size *yd)
2926{
2928 cpl_ensure_code (ia, CPL_ERROR_NULL_INPUT);
2929 cpl_ensure_code (ib, CPL_ERROR_NULL_INPUT);
2930 cpl_ensure_code (xd, CPL_ERROR_NULL_INPUT);
2931 cpl_ensure_code (yd, CPL_ERROR_NULL_INPUT);
2932
2933 cpl_error_code error = CPL_ERROR_NONE;
2934 cpl_type type = cpl_image_get_type(ia);
2935 cpl_size nx = cpl_image_get_size_x(ia);
2936 cpl_size ny = cpl_image_get_size_y(ia);
2937
2938 cpl_image * ic = cpl_image_new(nx, ny, type);
2939 cpl_image * fa = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2940 cpl_image * fb = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2941 cpl_image * fc = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2942
2943 cpl_imagelist * iab = cpl_imagelist_new();
2944 cpl_imagelist * fab = cpl_imagelist_new();
2945
2946 /* Put image ia into imagelist iab */
2947 error = cpl_imagelist_set(iab, ia, 0);
2948
2949 /* Put image ib into imagelist iab */
2950 error = cpl_imagelist_set(iab, ib, 1);
2951
2952 /* Put empty Fourier image fa in imagelist fab */
2953 error = cpl_imagelist_set(fab, fa, 0);
2954
2955 /* Put empty Fourier image fb in imagelist fab*/
2956 error = cpl_imagelist_set(fab, fb, 1);
2957
2958 /* FFT on the image list iab and fill FFT image list fab */
2959 error = cpl_fft_imagelist(fab, iab, CPL_FFT_FORWARD);
2960
2961 /* conjugate fourier image list fb to fc */
2962 error = cpl_image_conjugate(fc, fb);
2963
2964 /* Multiply fourier image list fa and fc*/
2965 error = cpl_image_multiply(fc, fa);
2966
2967 /* Reverse FFT */
2968 error = cpl_fft_image(ic, fc, CPL_FFT_BACKWARD | CPL_FFT_NOSCALE);
2969
2970 /* Find max position -> shift*/
2971 error = cpl_image_get_maxpos(ic, xd, yd);
2972
2973 /* Unwrap for "negative" maximum position */
2974 if (*xd > nx/2) *xd = *xd - nx;
2975 if (*yd > ny/2) *yd = *yd - ny;
2976
2977 /* subtract one for come from fits convention (lower left pixel = 1,1) to shift vector */
2978 (*xd)--;
2979 (*yd)--;
2980
2981 // cpl_image_save (ia, "ia.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2982 // cpl_image_save (ib, "ib.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2983 // cpl_image_save (ic, "ic.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2984
2985 /* Free memory allocated by this routine */
2986 FREE (cpl_imagelist_unwrap, iab);
2987 FREE (cpl_imagelist_unwrap, fab);
2988 FREE (cpl_image_delete, ic);
2989 FREE (cpl_image_delete, fa);
2990 FREE (cpl_image_delete, fb);
2991 FREE (cpl_image_delete, fc);
2992
2994 return (error);
2995}
2996
#define GRAVI_PRIMARY_HDR_EXT
Definition: gravi-test.h:212
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
Definition: gravi_data.h:39
#define gravi_data_get_header(data)
Definition: gravi_data.h:75
#define gravi_data_get_img(data, ext)
Definition: gravi_data.h:79
const cpl_size ntel
cpl_image_multiply(flat_profiled, profile_crop)
cpl_msg_debug(cpl_func, "Spectra has <50 pixels -> don't flat")
cpl_propertylist * header
Definition: gravi_old.c:2004
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
cpl_propertylist_update_double(header, "ESO QC MINWAVE SC", cpl_propertylist_get_double(plist, "ESO QC MINWAVE SC"))
cpl_image_delete(flat_profiled)
#define GRAVI_OI_VIS_ACQ_EXT
Definition: gravi_pfits.h:89
#define SINGLE_STS
Definition: gravi_pfits.h:182
#define DUAL_STS
Definition: gravi_pfits.h:183
#define INSNAME_ACQ
Definition: gravi_pfits.h:199
#define GRAVI_IMAGING_DATA_ACQ_EXT
Definition: gravi_pfits.h:41
#define MODE_SINGLE
Definition: gravi_pfits.h:177
#define MODE_ONAXIS
Definition: gravi_pfits.h:179
#define gravi_msg_function_exit(flag)
Definition: gravi_utils.h:85
#define FREE(function, variable)
Definition: gravi_utils.h:69
#define gravi_msg_function_start(flag)
Definition: gravi_utils.h:84
#define CPLCHECK_MSG(msg)
Definition: gravi_utils.h:45
#define GRAVI_MATH_RAD_MAS
Definition: gravi_utils.h:119
#define CPLCHECK(msg)
Definition: gravi_utils.h:42
cpl_error_code gravi_image_fft_correlate(cpl_image *ia, cpl_image *ib, cpl_size *xd, cpl_size *yd)
Correlate two images using FFT.
cpl_image * gravi_image_extract(cpl_image *image_in, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
extract sub window of image (similar to cpl_image_extract)
double exp1(double x)
Definition: gravi_acqcam.c:291
#define GRAVI_SPOT_NSPOT
Definition: gravi_acqcam.c:84
#define GRAVI_SPOT_NSEARCH
Definition: gravi_acqcam.c:86
#define GRAVI_SPOT_FWHM
Definition: gravi_acqcam.c:94
int gravi_acqcam_xy_diode(const double v[], double *xd, double *yd)
Definition: gravi_acqcam.c:301
cpl_error_code gravi_acqcam_set_pupil_table_v2(cpl_table *acqcam_table, cpl_propertylist *header, cpl_propertylist *o_header, cpl_array *good_frames, cpl_vector *scale_vector, cpl_array *bad_frames_short, cpl_bivector **diode_pos_offset, cpl_vector *focus_value, gravi_data *static_param_data)
store the pupil offset in the table, and create QC parameters.
#define GRAVI_SPOT_SWINDOW
Definition: gravi_acqcam.c:87
cpl_error_code gravi_acq_measure_max(cpl_image *img, double x, double y, double size, double *img_max)
Definition: gravi_acqcam.c:485
cpl_error_code gravi_acq_fit_gaussian(cpl_image *img, double *x, double *y, double *ex, double *ey, cpl_size size)
Fit a Gaussian into an image, and mark the position.
Definition: gravi_acqcam.c:516
#define GRAVI_SPOT_NFOCUS
Definition: gravi_acqcam.c:85
cpl_error_code gravi_acqcam_field(cpl_image *mean_img, cpl_imagelist *acqcam_imglist, cpl_propertylist *header, cpl_table *acqcam_table, cpl_propertylist *o_header, gravi_data *static_param_data)
Reduce the images of field from ACQ.
cpl_error_code gravi_acqcam_get_diode_ref_v2(cpl_propertylist *header, cpl_array *good_frames, cpl_vector *scale_vector, cpl_bivector **diode_pos_telescope, int nrow_on)
Get the position of the telescope diodes in pixels.
#define CHECK_BIT(var, pos)
Definition: gravi_acqcam.c:99
const int * GRAVI_LVMQ_FREE
Definition: gravi_acqcam.c:185
cpl_error_code gravi_acqcam_select_good_frames_v2(cpl_imagelist *acqcam_imglist, cpl_imagelist *pupilImage_onFrames, cpl_array *good_frames)
select pupil frames with pupil beacon on. Clean pupil frames by substraction of images with pupil bea...
Definition: gravi_acqcam.c:947
cpl_error_code gravi_acqcam_get_pup_ref_v2(cpl_propertylist *header, cpl_bivector *diode_pos_subwindow)
Get the reference pixels for the pupil guiding on the acquisition camera.
#define GRAVI_SPOT_NTEL
Definition: gravi_acqcam.c:82
#define GRAVI_SPOT_DIODE
Definition: gravi_acqcam.c:93
double sin1(double x)
Definition: gravi_acqcam.c:276
cpl_error_code gravi_acqcam_spot_imprint_v2(cpl_image *mean_img, cpl_bivector **diode_pos_offset, cpl_bivector **diode_pos_theoretical, int ury)
imprint zeros at the postion of the detected pupil spots
#define GRAVI_SPOT_NLENS
Definition: gravi_acqcam.c:83
cpl_error_code gravi_acqcam_perform_shiftandadd_v2(cpl_imagelist *pupilImage_onFrames, cpl_imagelist **pupilImage_shiftandadd, cpl_array *good_frames, cpl_vector *focus_value, cpl_bivector **diode_pos_theoretical, cpl_bivector **diode_pos_offset, cpl_size nrow_on)
perform shift and add of the pupil image according to the theoretical positions
#define GRAVI_SPOT_NA
Definition: gravi_acqcam.c:89
#define GRAVI_SPOT_FLUX
Definition: gravi_acqcam.c:95
double gravi_acqcam_z2meter(double PositionPixels, gravi_data *static_param_data)
Convert z_shift from [pixel] to [meters] Formula extracted from gvacqPupilTracker....
#define GRAVI_SPOT_ANGLE
Definition: gravi_acqcam.c:91
cpl_error_code gravi_acqcam_get_diode_theoretical_v2(cpl_bivector *diode_pos_subwindow, cpl_bivector **diode_pos_telescope, cpl_bivector **diode_pos_theoretical, cpl_size nrow_on, int ury)
Get the position of the telescope diodes in pixels.
cpl_error_code gravi_preproc_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *bad_map)
Preprocess the ACQ images: correct bad pixels, clean from pupil background via blinking,...
Definition: gravi_acqcam.c:208
cpl_error_code gravi_acqcam_clean_pupil_v2(cpl_imagelist *acqcam_imglist, cpl_imagelist *pupilImage_filtered, const cpl_size ury)
Cleaning pupil images by cross-correlation with gaussian function.
Definition: gravi_acqcam.c:818
cpl_error_code gravi_acq_measure_strehl(cpl_image *img, double x, double y, double pscale, double *SR, cpl_propertylist *header)
measure Strehl Ratio of the source at the given location
Definition: gravi_acqcam.c:423
cpl_error_code gravi_acqcam_pupil(cpl_image *mean_img, cpl_imagelist *acqcam_imglist, cpl_propertylist *header, cpl_table *acqcam_table, cpl_propertylist *o_header, gravi_data *static_param_data)
cpl_error_code gravi_acqcam_get_pupil_offset_v2(cpl_imagelist **pupilImage_shiftandadd, cpl_array *good_frames, cpl_bivector **diode_pos_offset, cpl_propertylist *o_header, cpl_size nrow_on)
on the images shift and added, find maximum and perform gaussian fit
#define GRAVI_SPOT_SCALE
Definition: gravi_acqcam.c:92
double gravi_acqcam_defocus_scaling(int focus)
gives focus value as in a list of value
cpl_error_code gravi_acqcam_pupil_v2(cpl_image *mean_img, cpl_imagelist *acqcam_imglist, cpl_propertylist *header, cpl_table *acqcam_table, cpl_propertylist *o_header, gravi_data *static_param_data)
Reduce the images of pupil from ACQ V2.0.
Definition: gravi_acqcam.c:634
#define GRAVI_SPOT_SUB
Definition: gravi_acqcam.c:90
cpl_error_code gravi_reduce_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *sky_data, gravi_data *dark_data, gravi_data *static_param_data, int saveAcqTable)
Reduce the ACQ camera images.
cpl_error_code gravi_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
Definition: gravi_cpl.c:1655
cpl_error_code gravi_image_replace_window(cpl_image *img1, const cpl_image *img2, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury, cpl_size llx2, cpl_size lly2)
Definition: gravi_cpl.c:2003
double gravi_table_get_column_std(cpl_table *table, const char *name, int base, int nbase)
Definition: gravi_cpl.c:387
cpl_error_code gravi_data_add_cube(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_imagelist *imglist)
Add an IMAGE (imagelist) extension in gravi_data.
Definition: gravi_data.c:2353
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
Definition: gravi_data.c:2049
cpl_error_code gravi_data_add_table(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_table *table)
Add a BINTABLE extension in gravi_data.
Definition: gravi_data.c:2289
cpl_error_code gravi_data_add_img(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_image *image)
Add an IMAGE (single image) extension in gravi_data.
Definition: gravi_data.c:2394
cpl_imagelist * gravi_data_get_cube(gravi_data *self, const char *extname)
Return a pointer on an IMAGE extension by its EXTNAME.
Definition: gravi_data.c:2131
int gravi_data_has_extension(gravi_data *raw_calib, const char *ext_name)
Check if data has extension with given EXTNAME.
Definition: gravi_data.c:1808
cpl_table * gravi_data_get_table(gravi_data *self, const char *extname)
Return a pointer on a table extension by its EXTNAME.
Definition: gravi_data.c:2096
double gravi_pfits_get_zenithangle_beamb_acqcam(const cpl_propertylist *plist, int tel, int n, int nrow)
Definition: gravi_pfits.c:850
double gravi_pfits_get_sobj_y(const cpl_propertylist *plist)
Definition: gravi_pfits.c:461
int gravi_pfits_get_mode(const cpl_propertylist *plist)
Definition: gravi_pfits.c:193
double gravi_pfits_get_northangle_acqcam(const cpl_propertylist *plist, int tel)
Definition: gravi_pfits.c:794
int gravi_pfits_get_axis(const cpl_propertylist *plist)
Definition: gravi_pfits.c:213
const char * gravi_pfits_get_feed(const cpl_propertylist *plist)
Definition: gravi_pfits.c:435
double gravi_pfits_get_sobj_x(const cpl_propertylist *plist)
Definition: gravi_pfits.c:455
double gravi_pfits_get_time_acqcam(const cpl_propertylist *header, cpl_size row)
Time of the middle of the ACQCAM exposure row in [us], counted from PRC.ACQ.START.
Definition: gravi_pfits.c:749
double gravi_pfits_get_ptfc_acqcam(const cpl_propertylist *plist, int spot)
Definition: gravi_pfits.c:773
cpl_error_code gravi_msg_warning(const char *component, const char *msg)
Definition: gravi_utils.c:127
const char * gravi_conf_get_telname(int gravi_beam, cpl_propertylist *header)
Definition: gravi_utils.c:1191