GRAVI Pipeline Reference Manual 1.9.2
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{
639cpl_ensure_code(acqcam_imglist, CPL_ERROR_NULL_INPUT);
640cpl_ensure_code(header, CPL_ERROR_NULL_INPUT);
641cpl_ensure_code(acqcam_table, CPL_ERROR_NULL_INPUT);
642cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
643
644cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (acqcam_imglist, 0));
645cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (acqcam_imglist, 0));
646
647/*
648 * constant ury is the beginning of the upper part of the acquisition camera
649 * It delimitates where the pupil beacon images starts
650 * Not sure it works in full frame case (ny>1100)
651 */
652const cpl_size ury = (ny>1100) ? 1200 : 745;
653
654/*
655 * First step. The goal is to get the images pupil, cutted to keep only pupil frames (using ury).
656 * The output is the pupil, after kernel filtering with a Gaussian pdf the size of the pupil beacon
657 */
658 cpl_imagelist * pupilImage_filtered = cpl_imagelist_new ();
659 gravi_acqcam_clean_pupil_v2(acqcam_imglist, pupilImage_filtered, ury);
660 CPLCHECK_MSG("Cannot clean pupil");
661
662/*
663 * Second step. The goal is to select the frames with the blinking beacons on
664 * The same routines gets the frames where the beacons are off
665 * It uses theses frames to remove the background
666 */
667 cpl_array * good_frames= cpl_array_new(cpl_imagelist_get_size(pupilImage_filtered),CPL_TYPE_INT);
668 cpl_imagelist * pupilImage_onFrames = cpl_imagelist_new ();
669 gravi_acqcam_select_good_frames_v2(pupilImage_filtered,pupilImage_onFrames,good_frames);
670 CPLCHECK_MSG("Cannot find blinking pupil frames");
671
672/* Third step. The goal is to initialize the vectors which will be used to know where the diodes
673 * are on the detectors (hence the use of bivectors, which contains x and y positions
674 * the size of the arrays are defined using several fixed parameters:
675 * GRAVI_SPOT_NTEL = 4 corresponds to the 4 telescopes
676 * GRAVI_SPOT_NSPOT = 4 corresponds to the 4 diodes on each telescopes
677 * GRAVI_SPOT_NLENS = 4 corresponds to the 4 lenslets on the acquisition camera
678 * nrow_on= 122 (for example) which corresponds to the N frames where the diodes are ON */
679
680 cpl_size nrow = cpl_imagelist_get_size(pupilImage_filtered);
681 cpl_size nrow_on = cpl_imagelist_get_size(pupilImage_onFrames);
682
683 if ((nrow_on<1)||(nrow<4))
684 {
685 /* exiting with empty pupil table */
686 cpl_array_fill_window_int (good_frames, 0, nrow_on, 0);
687 gravi_acqcam_set_pupil_table_v2(acqcam_table, header, o_header, good_frames, NULL, NULL, NULL, NULL, NULL);
688 CPLCHECK_MSG("Cannot stor pupil offset values in OI_ACQ table");
689
690 cpl_imagelist_delete(pupilImage_onFrames);
691 cpl_imagelist_delete(pupilImage_filtered);
692 cpl_array_delete(good_frames);
693 cpl_msg_warning (cpl_func, "Cannot reduce the %lli frames (not enough frames)",nrow);
695 return CPL_ERROR_NONE;
696 }
697
698 cpl_vector * scale_vector = cpl_vector_new(GRAVI_SPOT_NTEL);
699 cpl_vector * focus_value = cpl_vector_new(GRAVI_SPOT_NTEL);
700 cpl_array * bad_frames_short= cpl_array_new(nrow_on,CPL_TYPE_INT);
701 cpl_array_fill_window_int (bad_frames_short, 0, nrow_on, 0);
702 cpl_bivector * diode_pos_subwindow = cpl_bivector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NTEL);
703 cpl_bivector ** diode_pos_telescope = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
704 cpl_bivector ** diode_pos_theoretical = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
705 cpl_bivector ** diode_pos_offset = cpl_malloc (nrow_on * sizeof(cpl_bivector *)) ;
706 cpl_imagelist ** pupilImage_shiftandadd = cpl_malloc (GRAVI_SPOT_NTEL * sizeof(cpl_imagelist *));
707
708 for (int n = 0 ; n < nrow_on; n++)
709 {
710 diode_pos_telescope[n] = cpl_bivector_new (GRAVI_SPOT_NSPOT*GRAVI_SPOT_NTEL);
711 diode_pos_theoretical[n] = cpl_bivector_new (GRAVI_SPOT_NFOCUS*GRAVI_SPOT_NSPOT*GRAVI_SPOT_NLENS*GRAVI_SPOT_NTEL);
712 diode_pos_offset[n] = cpl_bivector_new (GRAVI_SPOT_NTEL);
713 }
714
715 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
716 pupilImage_shiftandadd[tel] = cpl_imagelist_new();
717 CPLCHECK_MSG("Cannot initialized vectors");
718
719/*
720 * Fourth step. Here we get the position of the 4 subwindows (correspond to each lenslet), in pixel/detector space
721 */
722 gravi_acqcam_get_pup_ref_v2(header, diode_pos_subwindow);
723 CPLCHECK_MSG("Cannot get pupil reference values");
724
725/*
726 * Fifth step. Here we get the position of the 4 diodes of each 4 telescopes, in pixel/detector space
727 */
728 gravi_acqcam_get_diode_ref_v2(header, good_frames, scale_vector, diode_pos_telescope, nrow_on);
729 CPLCHECK_MSG("Cannot get telescope beacons position");
730
731/*
732 * Seventh step. Get a large array of diode vector position (diode_pos_theoretical)
733 * this for each telescope, diode, lenslet, and focal number.
734 */
735 gravi_acqcam_get_diode_theoretical_v2(diode_pos_subwindow ,diode_pos_telescope,
736 diode_pos_theoretical, nrow_on, ury);
737 CPLCHECK_MSG("Cannot get theoretical position of beacons");
738
739/*
740 * Eigth step. Test each focus position to see which one is best. Then, collapse all spot images according to that focus value.
741 * The result pupilImage_shiftandadd correspond to an image for each blinking frame and each telescope
742 */
743 gravi_acqcam_perform_shiftandadd_v2(pupilImage_onFrames, pupilImage_shiftandadd, good_frames, focus_value, diode_pos_theoretical,
744 diode_pos_offset, nrow_on);
745 CPLCHECK_MSG("Cannot find correct focus to shift and add pupil images");
746
747/*
748 * Nineth step. For each one of the pupilImage_shiftandadd, fit a gaussian to know where the pupil is
749 * The measurement is done as an offset with respect to the theroetical position of the beacons
750 */
751 gravi_acqcam_get_pupil_offset_v2(pupilImage_shiftandadd, bad_frames_short, diode_pos_offset , o_header, nrow_on);
752 CPLCHECK_MSG("Cannot get pupil offset values");
753
754/*
755 * Last step. Storing the value of the pupil position in a table
756 */
757 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);
758 CPLCHECK_MSG("Cannot stor pupil offset values in OI_ACQ table");
759
760 /*
761 * Bonus. Add collapsed image to mean image and imprint best position as a cross in image
762 */
763
764 cpl_image * pupilImage_onFrames_collapse=cpl_imagelist_collapse_create (pupilImage_onFrames);
765 cpl_image * pupilImage_onFrames_1frame=cpl_imagelist_get (pupilImage_onFrames,0);
766 gravi_image_replace_window (mean_img, pupilImage_onFrames_1frame, 1, ury, nx, ny, 1, 1);
767 gravi_acqcam_spot_imprint_v2(mean_img, diode_pos_offset, diode_pos_theoretical, ury);
768 CPLCHECK_MSG("Cannot modify mean_image");
769
770/* cleaning all arrays ( free memory ) */
771
772for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
773{
774 if (pupilImage_shiftandadd[tel] != NULL)
775 cpl_imagelist_delete( pupilImage_shiftandadd[tel]);
776}
777FREE (cpl_free, pupilImage_shiftandadd);
778
779
780 for (int n = 0 ; n < nrow_on; n++)
781 {
782 cpl_bivector_delete( diode_pos_telescope[n]);
783 cpl_bivector_delete( diode_pos_theoretical[n]);
784 cpl_bivector_delete( diode_pos_offset[n]);
785 }
786 FREE (cpl_free, diode_pos_telescope);
787 FREE (cpl_free, diode_pos_theoretical);
788 FREE (cpl_free, diode_pos_offset);
789
790 cpl_bivector_delete( diode_pos_subwindow);
791 cpl_image_delete(pupilImage_onFrames_collapse);
792 cpl_imagelist_delete(pupilImage_onFrames);
793 cpl_imagelist_delete(pupilImage_filtered);
794 cpl_array_delete(good_frames);
795 cpl_array_delete(bad_frames_short);
796 cpl_vector_delete(scale_vector);
797 cpl_vector_delete(focus_value);
798
799CPLCHECK_MSG("Failed at freeing memory");
801return CPL_ERROR_NONE;
802}
803
804/*----------------------------------------------------------------------------*/
815/*----------------------------------------------------------------------------*/
816
817cpl_error_code gravi_acqcam_clean_pupil_v2(cpl_imagelist * acqcam_imglist, cpl_imagelist * pupilImage_filtered, const cpl_size ury)
818{
820
821 cpl_ensure_code(acqcam_imglist, CPL_ERROR_NULL_INPUT);
822 cpl_ensure_code(pupilImage_filtered, CPL_ERROR_NULL_INPUT);
823 cpl_ensure_code(ury, CPL_ERROR_NULL_INPUT);
824
825 /* Number of row */
826 cpl_size nrow = cpl_imagelist_get_size(acqcam_imglist);
827
828 cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (acqcam_imglist, 0));
829 cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (acqcam_imglist, 0));
830
831 cpl_imagelist * pupilImage = cpl_imagelist_new ();
832
833 for (cpl_size n = 0; n < nrow ; n++)
834 {
835 cpl_image * small_img_tmp = cpl_imagelist_get (acqcam_imglist, n);
836 cpl_image * small_img = cpl_image_extract( small_img_tmp, 1, ury, nx, ny);
837 cpl_imagelist_set(pupilImage, small_img,n);
838 }
839
840 nx=cpl_image_get_size_x (cpl_imagelist_get (pupilImage, 0));
841 ny=cpl_image_get_size_y (cpl_imagelist_get (pupilImage, 0));
842
843 cpl_image* pupilImageBackground = cpl_imagelist_collapse_minmax_create (pupilImage,(cpl_size)(nrow/8),(cpl_size)(5*nrow/8));
844
845 cpl_imagelist_subtract_image (pupilImage,pupilImageBackground);
846 CPLCHECK_MSG("Failure to subtract background");
847
848 cpl_size Npref = 21;
849 cpl_matrix * kernel1 = cpl_matrix_new (Npref,Npref);
850 cpl_matrix * kernel2 = cpl_matrix_new (Npref,Npref);
851 for (cpl_size x =0;x< Npref; x++)
852 for (cpl_size y =0;y< Npref; y++)
853 {
854 double radius_square=(x-(Npref-1)/2)*(x-(Npref-1)/2)+(y-(Npref-1)/2)*(y-(Npref-1)/2);
855 cpl_matrix_set(kernel1,x,y,exp(- radius_square /30));
856 cpl_matrix_set(kernel2,x,y,exp(- radius_square /40));
857 }
858
859 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_mean (kernel1));
860 cpl_matrix_divide_scalar (kernel2, cpl_matrix_get_mean (kernel2));
861 cpl_matrix_subtract (kernel1,kernel2);
862 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_stdev (kernel1));
863 cpl_matrix_divide_scalar (kernel1, cpl_matrix_get_max (kernel1));
864
865 CPLCHECK_MSG("Error computing gaussian Kernel");
866
867 // zero-pad to match image and kernel size, plus extra buffer for spurious data at edges
868 cpl_size fft_nx = nx + Npref - 1;
869 cpl_size fft_ny = ny + Npref - 1;
870 cpl_size half_kernel_w = (Npref - 1) / 2;
871
872 // Copy kernel into image with zero-padding
873 cpl_matrix_resize(kernel1, 0, ny-1, 0, nx-1);
874 cpl_image *kernel_image = cpl_image_wrap_double(fft_nx, fft_ny, cpl_matrix_get_data(kernel1));
875
876 // Compute kernel FFT
877 cpl_image *kernel_fft = cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE_COMPLEX);
878 cpl_fft_image(kernel_fft, kernel_image, CPL_FFT_FORWARD);
879
880 // Create imagelist with zero-padded copies of the source images
881 cpl_imagelist *pupilImage_padded = cpl_imagelist_new();
882 for (cpl_size n = 0; n < nrow ; n++)
883 {
884 const cpl_image *small_img_tmp = cpl_imagelist_get_const(pupilImage, n);
885 cpl_image *padded_img_tmp = cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE);
886 cpl_image_copy(padded_img_tmp, small_img_tmp, 1, 1);
887 cpl_imagelist_set(pupilImage_padded, padded_img_tmp, n);
888 }
889
890 // Fourier-transform the images
891 cpl_imagelist *pupilImage_fft = cpl_imagelist_new();
892 for (cpl_size n = 0; n < nrow; n++)
893 cpl_imagelist_set(pupilImage_fft, cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE_COMPLEX), n);
894 cpl_fft_imagelist(pupilImage_fft, pupilImage_padded, CPL_FFT_FORWARD | CPL_FFT_FIND_MEASURE);
895
896 // Convolve with kernel
897 cpl_imagelist_multiply_image(pupilImage_fft, kernel_fft);
898
899 // Reverse FFT to get filtered images
900 cpl_imagelist *pupilImage_filtered_padded = cpl_imagelist_new();
901 for (cpl_size n = 0; n < nrow; n++)
902 cpl_imagelist_set(pupilImage_filtered_padded, cpl_image_new(fft_nx, fft_ny, CPL_TYPE_DOUBLE), n);
903 cpl_fft_imagelist(pupilImage_filtered_padded, pupilImage_fft, CPL_FFT_BACKWARD | CPL_FFT_FIND_MEASURE);
904
905 // Trim the padding from the image edges
906 for (cpl_size n = 0; n < nrow; n++) {
907 const cpl_image *pad_filt_img_tmp = cpl_imagelist_get_const(pupilImage_filtered_padded, n);
908 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);
909 cpl_imagelist_set(pupilImage_filtered, trimmed_img, n);
910 }
911
912 cpl_matrix_delete(kernel1);
913 cpl_matrix_delete(kernel2);
914
915 cpl_image_unwrap(kernel_image);
916 cpl_image_delete(kernel_fft);
917
918 cpl_image_delete(pupilImageBackground);
919 cpl_imagelist_delete(pupilImage);
920 cpl_imagelist_delete(pupilImage_fft);
921 cpl_imagelist_delete(pupilImage_padded);
922 cpl_imagelist_delete(pupilImage_filtered_padded);
923
924 CPLCHECK_MSG("Pupil Fitting V2 does not work");
926 return CPL_ERROR_NONE;
927}
928
929
930/*----------------------------------------------------------------------------*/
944/*----------------------------------------------------------------------------*/
945
946cpl_error_code gravi_acqcam_select_good_frames_v2(cpl_imagelist * pupilImage_filtered, cpl_imagelist * pupilImage_onFrames, cpl_array * good_frames)
947{
949 int nv =0;
950
951 cpl_ensure_code(pupilImage_filtered, CPL_ERROR_NULL_INPUT);
952 cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
953
954 /* Number of row , size of image */
955 cpl_size nrow = cpl_imagelist_get_size(pupilImage_filtered);
956 cpl_size nx = cpl_image_get_size_x (cpl_imagelist_get (pupilImage_filtered, 0));
957 cpl_size ny = cpl_image_get_size_y (cpl_imagelist_get (pupilImage_filtered, 0));
958
959 cpl_vector * pupil_max = cpl_vector_new(4);
960 cpl_vector * pupil_mediam = cpl_vector_new(nrow);
961 cpl_array * back_frames= cpl_array_new(nrow,CPL_TYPE_INT);
962
963 for (cpl_size n = 0; n < nrow; n++)
964 {
965
966 cpl_image * image = cpl_imagelist_get (pupilImage_filtered,n);
967 cpl_vector_set(pupil_max, 0, cpl_image_get_max_window (image, 1, 1, 250, ny));
968 cpl_vector_set(pupil_max, 1, cpl_image_get_max_window (image, 251, 1, 500, ny));
969 cpl_vector_set(pupil_max, 2, cpl_image_get_max_window (image, 501, 1, 750, ny));
970 cpl_vector_set(pupil_max, 3, cpl_image_get_max_window (image, 751, 1, nx, ny));
971 cpl_vector_set(pupil_mediam, n, cpl_vector_get_median (pupil_max));
972 /*
973 cpl_msg_info (cpl_func, "tototo image valeur %lli = %.2f",n, cpl_image_get(image, 101, 151, &nv));
974 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));
975 cpl_msg_info (cpl_func, "Finding min brightness level for beacon ON %lli: %.5f ADU",n, cpl_vector_get(pupil_mediam, n));*/
976
977 }
978
979 cpl_vector * pupil_mediam_sort = cpl_vector_duplicate (pupil_mediam);
980 cpl_vector_sort (pupil_mediam_sort,CPL_SORT_ASCENDING);
981 cpl_vector * pupil_median_down = cpl_vector_extract (pupil_mediam_sort, 0, (nrow-1)/2, 1);
982 cpl_vector * pupil_median_up = cpl_vector_extract (pupil_mediam_sort, nrow/2, nrow-1, 1);
983 double threshold_up = cpl_vector_get_median(pupil_median_up);
984 double threshold_down = cpl_vector_get_median(pupil_median_down);
985 double threshold=(threshold_up+threshold_down)/2.0;
986
987 CPLCHECK_MSG("Could not get median maximum of beacon flux ");
988 cpl_msg_info (cpl_func, "Found threshold brightness level for beacon : %.5f ADU", threshold);
989
990 for (cpl_size n = 0; n < nrow; n++)
991 {
992 if (cpl_vector_get(pupil_mediam, n) >threshold *1.1)
993 cpl_array_set_int (good_frames,n,1);
994 else
995 cpl_array_set_int (good_frames,n,0);
996 if (cpl_vector_get(pupil_mediam, n) <threshold * 0.9)
997 cpl_array_set_int (back_frames,n,1);
998 else
999 cpl_array_set_int (back_frames,n,0);
1000 }
1001
1002 cpl_size frames_background = 4;
1003 cpl_size n_goodFrames = 0;
1004
1005 CPLCHECK_MSG("Failed to find blincking pupil files");
1006
1007 for (cpl_size n = 0; n < nrow; n++)
1008 if (cpl_array_get_int(good_frames,n,&nv) == 1)
1009 {
1010 cpl_image * image = cpl_imagelist_get (pupilImage_filtered,n);
1011 /*cpl_msg_info (cpl_func, "Good frame image valeur %lli = %.2f",n, cpl_image_get(image, 101, 151, &nv));*/
1012
1013 cpl_size n_frames_background = 0;
1014 cpl_imagelist * pupilImage_offFrames = cpl_imagelist_new ();
1015 for (cpl_size b = 1; b < nrow; b++)
1016 {
1017 if (n+b < nrow)
1018 if ((cpl_array_get_int(back_frames,n+b,&nv) == 1)&(n_frames_background<frames_background))
1019 {
1020 cpl_image * image_background = cpl_imagelist_get (pupilImage_filtered,n+b);
1021 cpl_image * image_background_copy = cpl_image_duplicate (image_background);
1022 cpl_imagelist_set(pupilImage_offFrames, image_background_copy, n_frames_background);
1023 n_frames_background = n_frames_background + 1;
1024 }
1025 if (n-b >= 0)
1026 if ((cpl_array_get_int(back_frames,n-b,&nv) == 1)&(n_frames_background<frames_background))
1027 {
1028 cpl_image * image_background = cpl_imagelist_get (pupilImage_filtered,n-b);
1029 cpl_image * image_background_copy = cpl_image_duplicate (image_background);
1030 cpl_imagelist_set(pupilImage_offFrames, image_background_copy, n_frames_background);
1031 n_frames_background = n_frames_background + 1;
1032 }
1033 }
1034
1035 cpl_image * image_background_mean = cpl_imagelist_collapse_create (pupilImage_offFrames);
1036 cpl_image * image_background_subtracted = cpl_image_subtract_create (image, image_background_mean);
1037 cpl_imagelist_set(pupilImage_onFrames, image_background_subtracted, n_goodFrames);
1038
1039 cpl_image_delete(image_background_mean);
1040 cpl_imagelist_delete(pupilImage_offFrames);
1041
1042 n_goodFrames ++;
1043 }
1044 cpl_msg_info(cpl_func, "Found %lli frames with pupil beacon ON (over %lli frames)",n_goodFrames, nrow);
1045
1046 cpl_vector_delete(pupil_max);
1047 cpl_vector_delete(pupil_mediam);
1048 cpl_vector_delete(pupil_median_up);
1049 cpl_vector_delete(pupil_median_down);
1050 cpl_vector_delete(pupil_mediam_sort);
1051 cpl_array_delete(back_frames);
1052
1053 CPLCHECK_MSG("Pupil Selection of good files failed");
1055 return CPL_ERROR_NONE;
1056}
1057
1058/*----------------------------------------------------------------------------*/
1070/*----------------------------------------------------------------------------*/
1071
1072cpl_error_code gravi_acqcam_get_pup_ref_v2 (cpl_propertylist * header, cpl_bivector * diode_pos_subwindow)
1073{
1075 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1076 cpl_ensure_code (diode_pos_subwindow, CPL_ERROR_NULL_INPUT);
1077
1078 cpl_size nsx = 0, nsy = 0, sx = 0, sy = 0;
1079
1080 cpl_vector * x_pos_subwindow = cpl_bivector_get_x (diode_pos_subwindow);
1081 cpl_vector * y_pos_subwindow = cpl_bivector_get_y (diode_pos_subwindow);
1082
1083
1084 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1085 {
1086 /* If sub-windowing, we read the sub-window size and
1087 * the sub-window start for pupil */
1088 if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
1089 char name[90];
1090
1091 nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
1092 sprintf (name, "ESO DET1 FRAM%d STRX", 3*GRAVI_SPOT_NTEL + tel + 1);
1093 sx = cpl_propertylist_get_int (header, name);
1094
1095 nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY");
1096 sprintf (name, "ESO DET1 FRAM%d STRY", 3*GRAVI_SPOT_NTEL + tel + 1);
1097 sy = cpl_propertylist_get_int (header, name);
1098
1099 CPLCHECK_MSG ("Cannot get sub-windowing parameters");
1100 }
1101
1102 cpl_msg_info (cpl_func,"sub-window pupil %i sx= %lld sy = %lld", tel, sx, sy);
1103
1104 /* Read the sub-apperture reference positions
1105 * Converted to accound for sub-windowing
1106 * In vector convention (start at 0,0) */
1107
1108 for (int lens = 0; lens < GRAVI_SPOT_NLENS ; lens++) {
1109 double xv = gravi_pfits_get_ptfc_acqcam (header, lens*GRAVI_SPOT_NTEL + tel + 1);
1110 double yv = gravi_pfits_get_ptfc_acqcam (header, lens*GRAVI_SPOT_NTEL + tel + 17);
1111 cpl_vector_set( x_pos_subwindow, lens*GRAVI_SPOT_NTEL+tel , xv - (sx - tel*nsx) - 1);
1112 cpl_vector_set( y_pos_subwindow, lens*GRAVI_SPOT_NTEL+tel , yv - (sy - 3*nsy) - 1);
1113 /* cpl_msg_info (cpl_func,"pupil %lli subC %i = %10.4f,%10.4f",
1114 tel, lens, x_diode_subwindow[lens][tel], y_diode_subwindow[lens][tel]);*/
1115 CPLCHECK_MSG ("Cannot get pupil reference position");
1116 }
1117 }
1118
1120 return CPL_ERROR_NONE;
1121}
1122
1123
1124/*----------------------------------------------------------------------------*/
1139/*----------------------------------------------------------------------------*/
1140
1141
1142cpl_error_code gravi_acqcam_get_diode_ref_v2 (cpl_propertylist * header,
1143 cpl_array * good_frames,
1144 cpl_vector * scale_vector,
1145 cpl_bivector ** diode_pos_telescope,
1146 int nrow_on)
1147{
1149 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1150 cpl_ensure_code (good_frames, CPL_ERROR_NULL_INPUT);
1151
1152 int nv = 0;
1153 double dx,dy,scale;
1154 double parang1 = cpl_propertylist_get_double(header, "ESO ISS PARANG START");
1155 double parang2 = cpl_propertylist_get_double(header, "ESO ISS PARANG END");
1156 CPLCHECK_MSG("Cannot determine parallactic angle");
1157
1158 /* check the feed mode */
1159 int sts_mode = 0;
1160 if ( !strcmp(gravi_pfits_get_feed (header), "DUAL_STS") ) {
1161 sts_mode = DUAL_STS;
1162 }
1163 else{
1164 sts_mode = SINGLE_STS;
1165 }
1166
1167 cpl_size nrow = cpl_array_get_size(good_frames);
1168 if (fabs(cpl_array_get_mean(good_frames)*nrow -nrow_on) > 1e-2)
1169 cpl_msg_error (cpl_func, "Ratio of blinking frames different %f from %i",cpl_array_get_mean(good_frames)*nrow,nrow_on);
1170
1171 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1172 {
1173
1174 /* Get the telescope name and ID */
1175 const char * telname = gravi_conf_get_telname (tel, header);
1176
1177 /* Check telescope name */
1178 if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
1179 cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
1180 CPLCHECK_MSG("Cannot get telescope names");
1181
1182 /* Hardcoded theoretical positions in mm */
1183
1184 /* If UTs or ATs, select scaling, rotation, and spacing */
1185 if (telname[0] == 'U') {
1186 // FE 2019-05-15 replaced with median measured for the whole 2017-2018
1187 // Galactic Center data set, which should be the best information at hand
1188 // cpl_vector_set (output, GRAVI_SPOT_SCALE, 16.225);
1189 if ( tel == 0 ) scale = 16.83;
1190 if ( tel == 1 ) scale = 17.42;
1191 if ( tel == 2 ) scale = 16.85;
1192 if ( tel == 3 ) scale = 17.41;
1193 // below information could be calculated from diode position in header
1194 // this would also give the 45 offset angle introduced above to calculate
1195 // the spot angle
1196 dx = 0.363;
1197 dy = 0.823;
1198 } else if (telname[0] == 'A') {
1199 // FE 2019-05-15 maybe we should also update the AT numbers
1200 scale = 73.0154;
1201 // EW, FE 2019-09-11: short and long side of AT beacons are
1202 // flipped when compared to UTs
1203 dx = 0.158;
1204 dy = 0.122;
1205 } else
1206 return cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1207 "Cannot get telescope name");
1208 CPLCHECK_MSG("Cannot get telescope scale");
1209 cpl_vector_set(scale_vector,tel,scale);
1210
1211 /* We see that fitting the spot angle fails if one spot is fainter than reflection from brightest spot
1212 resulting in wrong pupil_x,y, therefore we calculate angle from header information
1213 by design of telescope spiders. This model angle should be
1214 angle = north_angle + paralactic angle + 45
1215 Remark: checking the angles for the GC easter flare night, we get an offset of 45.75 degrees */
1216 double northangle = 0;
1217 double zenithangle = 0;
1218 if (sts_mode == SINGLE_STS) {
1219 northangle = gravi_pfits_get_northangle_acqcam(header, tel);
1220 CPLCHECK_MSG("Cannot determine field angle");
1221 }
1222
1223
1224 cpl_size n_on= 0;
1225 for (cpl_size n =0; n < nrow; n++)
1226 {
1227 /*cpl_msg_info (cpl_func, "Totos: %lli, %.2f",cpl_array_get_int(good_frames,n,&nv),cpl_array_get_mean(good_frames));*/
1228 if (cpl_array_get_int(good_frames,n,&nv) ==1)
1229 {
1230 double padif = parang2 - parang1;
1231 if (padif > 180)
1232 padif -= 360;
1233 if (padif < -180)
1234 padif += 360;
1235 double parang= parang1 + padif*n/(nrow-1);
1236
1237 // TODO: FE handle case of calibration unit //
1238 double angle = 0;
1239 if (sts_mode == DUAL_STS) {
1240 /* northangle = gravi_pfits_get_northangle_beamb_acqcam(header, tel, n, nrow); */
1241 zenithangle = gravi_pfits_get_zenithangle_beamb_acqcam(header, tel, n, nrow);
1242 angle = zenithangle - 45.;
1243 } else {
1244 angle = northangle + parang + 45.;
1245 }
1246
1247 if (angle < 0) angle += 180;
1248 if (angle > 180) angle -= 180;
1249 double cang = cos(angle * CPL_MATH_RAD_DEG) * scale;
1250 double sang = sin(angle * CPL_MATH_RAD_DEG) * scale;
1251
1252 /* Diode arrangement */
1253
1254 cpl_vector * x_pos_telescope = cpl_bivector_get_x (diode_pos_telescope[n_on]);
1255 cpl_vector * y_pos_telescope = cpl_bivector_get_y (diode_pos_telescope[n_on]);
1256
1257 cpl_vector_set(x_pos_telescope, 0*GRAVI_SPOT_NTEL+tel, -cang * dx + sang * dy);
1258 cpl_vector_set(x_pos_telescope, 1*GRAVI_SPOT_NTEL+tel, -cang * dx - sang * dy);
1259 cpl_vector_set(x_pos_telescope, 2*GRAVI_SPOT_NTEL+tel, cang * dx - sang * dy);
1260 cpl_vector_set(x_pos_telescope, 3*GRAVI_SPOT_NTEL+tel, cang * dx + sang * dy);
1261
1262 cpl_vector_set(y_pos_telescope, 0*GRAVI_SPOT_NTEL+tel, cang * dy + sang * dx);
1263 cpl_vector_set(y_pos_telescope, 1*GRAVI_SPOT_NTEL+tel, -cang * dy + sang * dx);
1264 cpl_vector_set(y_pos_telescope, 2*GRAVI_SPOT_NTEL+tel, -cang * dy - sang * dx);
1265 cpl_vector_set(y_pos_telescope, 3*GRAVI_SPOT_NTEL+tel, cang * dy - sang * dx);
1266
1267 n_on ++;
1268
1269 if (n_on == 10)
1270 {
1271 cpl_msg_info (cpl_func, "angle, parang %.2f, %.2f",angle,parang);
1272 cpl_msg_info (cpl_func, "dx and dy %.2f, %.2f",cang,sang);
1273 }
1274
1275
1276 }
1277 }
1278
1279
1280 CPLCHECK_MSG("Cannot determine diode position in telescope space");
1281 }
1282
1284 return CPL_ERROR_NONE;
1285}
1286
1287/*----------------------------------------------------------------------------*/
1305/*----------------------------------------------------------------------------*/
1306
1307cpl_error_code gravi_acqcam_get_diode_theoretical_v2(cpl_bivector * diode_pos_subwindow,
1308 cpl_bivector ** diode_pos_telescope,
1309 cpl_bivector ** diode_pos_theoretical,
1310 cpl_size nrow_on, int ury)
1311{
1312
1314 cpl_ensure_code (diode_pos_subwindow, CPL_ERROR_NULL_INPUT);
1315 cpl_ensure_code (diode_pos_telescope, CPL_ERROR_NULL_INPUT);
1316 cpl_ensure_code (diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1317
1318 cpl_vector * x_pos_subwindow = cpl_bivector_get_x (diode_pos_subwindow);
1319 cpl_vector * y_pos_subwindow = cpl_bivector_get_y (diode_pos_subwindow);
1320
1321 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1322 {
1323 double x_lenslet_mean=0.0;
1324 double y_lenslet_mean=0.0;
1325 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1326 {
1327 x_lenslet_mean+=cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)/GRAVI_SPOT_NLENS;
1328 y_lenslet_mean+=cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)/GRAVI_SPOT_NLENS;
1329 }
1330
1331 cpl_msg_info (cpl_func, "Reference pixel position for tel %d : X = %.2f, Y= %.2f", tel, x_lenslet_mean, y_lenslet_mean);
1332
1333 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1334 for (int focus = 0 ; focus < GRAVI_SPOT_NFOCUS; focus++)
1335 {
1336 double xlenslet=cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)+
1337 (cpl_vector_get(x_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)-x_lenslet_mean)*gravi_acqcam_defocus_scaling(focus);
1338 double ylenslet=cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)+
1339 (cpl_vector_get(y_pos_subwindow,lens*GRAVI_SPOT_NTEL+tel)-y_lenslet_mean)*gravi_acqcam_defocus_scaling(focus);
1340
1341 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1342 for (cpl_size n_on = 0 ; n_on < nrow_on; n_on++)
1343 {
1344
1345 cpl_vector * x_pos_telescope = cpl_bivector_get_x (diode_pos_telescope[n_on]);
1346 cpl_vector * y_pos_telescope = cpl_bivector_get_y (diode_pos_telescope[n_on]);
1347 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n_on]);
1348 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n_on]);
1349
1350 double x_theo = xlenslet + cpl_vector_get(x_pos_telescope,spot*GRAVI_SPOT_NTEL+tel);
1351 double y_theo = ylenslet + cpl_vector_get(y_pos_telescope,spot*GRAVI_SPOT_NTEL+tel)- ury +1;
1352
1353 cpl_vector_set(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel,x_theo);
1354 cpl_vector_set(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel,y_theo);
1355 }
1356 }
1357 }
1358
1359CPLCHECK_MSG("Cannot compute the theoretical coordinates of spots");
1361return CPL_ERROR_NONE;
1362
1363}
1364
1365/*----------------------------------------------------------------------------*/
1375/*----------------------------------------------------------------------------*/
1376
1377
1378cpl_error_code gravi_acqcam_spot_imprint_v2(cpl_image *mean_img, cpl_bivector ** diode_pos_offset, cpl_bivector ** diode_pos_theoretical, int ury)
1379{
1381
1382 cpl_ensure_code(mean_img, CPL_ERROR_NULL_INPUT);
1383 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1384 cpl_ensure_code(diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1385
1386 cpl_size nx = cpl_image_get_size_x (mean_img);
1387 cpl_size ny = cpl_image_get_size_y (mean_img);
1388
1389 /* imprint position for first image */
1390 int n_on = 0;
1391 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n_on]);
1392 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n_on]);
1393 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n_on]);
1394 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n_on]);
1395
1396 /* imprint position with no defocus */
1397 int focus = GRAVI_SPOT_NFOCUS/2;
1398
1399 /* Loop on diode and appertures */
1400 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1401 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1402 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1403 {
1404 double x_off = cpl_vector_get(x_pos_offset,tel);
1405 double y_off = cpl_vector_get(y_pos_offset,tel)+ury;
1406 double x_diode = cpl_vector_get(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1407 double y_diode = cpl_vector_get(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1408 cpl_size xf = roundl(x_diode + x_off) + 1;
1409 cpl_size yf = roundl(y_diode + y_off) + 1;
1410 if (xf < 2 || xf > nx-2 || yf < 2 || yf > ny-2) continue;
1411 cpl_image_set (mean_img, xf, yf, 0);
1412 cpl_image_set (mean_img, xf-1, yf, 0);
1413 cpl_image_set (mean_img, xf+1, yf, 0);
1414 cpl_image_set (mean_img, xf, yf+1, 0);
1415 cpl_image_set (mean_img, xf, yf-1, 0);
1416 CPLCHECK ("Cannot imprint cross in image");
1417 }
1418
1419 CPLCHECK ("Cannot imprint cross in image");
1421 return CPL_ERROR_NONE;
1422}
1423
1424
1425/*----------------------------------------------------------------------------*/
1456/*----------------------------------------------------------------------------*/
1457
1458
1459cpl_error_code gravi_acqcam_perform_shiftandadd_v2(cpl_imagelist * pupilImage_onFrames, cpl_imagelist ** pupilImage_shiftandadd, cpl_array * good_frames,
1460 cpl_vector * focus_value,
1461 cpl_bivector ** diode_pos_theoretical ,
1462 cpl_bivector ** diode_pos_offset, cpl_size nrow_on)
1463{
1465
1466 cpl_ensure_code(pupilImage_onFrames, CPL_ERROR_NULL_INPUT);
1467 cpl_ensure_code(pupilImage_shiftandadd, CPL_ERROR_NULL_INPUT);
1468 cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
1469 cpl_ensure_code(focus_value, CPL_ERROR_NULL_INPUT);
1470 cpl_ensure_code(diode_pos_theoretical, CPL_ERROR_NULL_INPUT);
1471 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1472 cpl_ensure_code(nrow_on, CPL_ERROR_NULL_INPUT);
1473
1474 if (cpl_imagelist_get_size(pupilImage_onFrames) != nrow_on)
1475 cpl_msg_error (cpl_func, "Problem with number of blinking frames");
1476
1477 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1478 {
1479 cpl_vector * focus_max = cpl_vector_new(GRAVI_SPOT_NFOCUS);
1480 for (cpl_size focus = 0 ; focus < GRAVI_SPOT_NFOCUS; focus++)
1481 {
1482 cpl_imagelist * pupilImage_shifted= cpl_imagelist_new();
1483 cpl_size n_images = 0;
1484 for (cpl_size n = 0 ; n < nrow_on; n++)
1485 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1486 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1487 {
1488 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n]);
1489 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n]);
1490
1491 double x_theo=cpl_vector_get(x_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1492 double y_theo=cpl_vector_get(y_pos_theoretical,((focus*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1493
1494 cpl_size x=(int) (x_theo + 0.5);
1495 cpl_size y=(int) (y_theo+ 0.5);
1496 cpl_image * verysmall_img_tmp = cpl_imagelist_get (pupilImage_onFrames, n);
1497 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);
1498 cpl_imagelist_set(pupilImage_shifted, verysmall_img, n_images);
1499 n_images ++;
1500 }
1501
1502 cpl_image * image_mean = cpl_imagelist_collapse_create (pupilImage_shifted);
1503 cpl_vector_set(focus_max, focus, cpl_image_get_max(image_mean));
1504 cpl_imagelist_delete(pupilImage_shifted);
1505 cpl_image_delete(image_mean);
1506 }
1507
1508 int focus_max_pos = cpl_vector_get_maxpos (focus_max);
1509 cpl_msg_info (cpl_func, "Focus value for telescope %d : F = %.2f pixels", tel, GRAVI_SPOT_SWINDOW * gravi_acqcam_defocus_scaling(focus_max_pos));
1510 cpl_vector_set(focus_value,tel,gravi_acqcam_defocus_scaling(focus_max_pos)*100.);
1511 CPLCHECK_MSG("Cannot find optimum focus position");
1512
1513 for (cpl_size n = 0 ; n < nrow_on; n++)
1514 {
1515 cpl_vector * x_pos_theoretical = cpl_bivector_get_x (diode_pos_theoretical[n]);
1516 cpl_vector * y_pos_theoretical = cpl_bivector_get_y (diode_pos_theoretical[n]);
1517 cpl_imagelist * pupilImage_shifted= cpl_imagelist_new();
1518 cpl_vector * offs_x = cpl_vector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NSPOT);
1519 cpl_vector * offs_y = cpl_vector_new (GRAVI_SPOT_NLENS*GRAVI_SPOT_NSPOT);
1520 cpl_size n_images = 0;
1521 for (int lens = 0 ; lens < GRAVI_SPOT_NLENS; lens++)
1522 for (int spot = 0 ; spot < GRAVI_SPOT_NSPOT; spot++)
1523 {
1524 double x_theo=cpl_vector_get(x_pos_theoretical,((focus_max_pos*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1525 double y_theo=cpl_vector_get(y_pos_theoretical,((focus_max_pos*GRAVI_SPOT_NSPOT+spot)*GRAVI_SPOT_NLENS+lens)*GRAVI_SPOT_NTEL+tel);
1526 int x= (int) (x_theo+ 0.5);
1527 int y= (int) (y_theo + 0.5);
1528 cpl_image * verysmall_img_tmp = cpl_imagelist_get (pupilImage_onFrames, n);
1529 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);
1530 /*if (n==10)
1531 {
1532 cpl_msg_info (cpl_func, "Cutting at position : X/Y = %lli/%lli, [%lli,%lli,%lli]", x,y,lens,spot,tel);
1533 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));
1534 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));
1535 }*/
1536 /*if (n == 10)
1537 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]);*/
1538 cpl_imagelist_set(pupilImage_shifted, verysmall_img, n_images);
1539 cpl_vector_set (offs_x, n_images, (x-x_theo));
1540 cpl_vector_set (offs_y, n_images, (y-y_theo));
1541
1542 n_images ++;
1543 }
1544
1545 cpl_bivector * offs = cpl_bivector_wrap_vectors (offs_x, offs_y);
1546
1547 double ppos_x, ppos_y;
1548 cpl_image** cpl_image_combined = cpl_geom_img_offset_saa (pupilImage_shifted,offs,CPL_KERNEL_DEFAULT,
1549 0,0,CPL_GEOM_INTERSECT,&ppos_x,&ppos_y);
1550 CPLCHECK_MSG("Cannot do fine shift and add");
1551
1552 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n]);
1553 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n]);
1554
1555 cpl_image * image_mean = cpl_image_extract (cpl_image_combined[0], 4, 4, GRAVI_SPOT_NSEARCH*2-3, GRAVI_SPOT_NSEARCH*2-3);
1556
1557 cpl_vector_set(x_pos_offset, tel, ppos_x+3);
1558 cpl_vector_set(y_pos_offset, tel, ppos_y+3);
1559
1560 cpl_imagelist_set(pupilImage_shiftandadd[tel], image_mean, n);
1561 CPLCHECK_MSG("Cannot add shift and added image to imagelist");
1562
1563 if (cpl_image_combined[0] != NULL) cpl_image_delete(cpl_image_combined[0]);
1564 if (cpl_image_combined[1] != NULL) cpl_image_delete(cpl_image_combined[1]);
1565 cpl_free(cpl_image_combined);
1566
1567 cpl_bivector_delete(offs);
1568 cpl_imagelist_delete(pupilImage_shifted);
1569 }
1570 cpl_vector_delete(focus_max);
1571 }
1572
1573 CPLCHECK_MSG("Fail at running the shift and add ESO algorithm");
1575
1576 return CPL_ERROR_NONE;
1577}
1578
1579
1580/*----------------------------------------------------------------------------*/
1598/*----------------------------------------------------------------------------*/
1599
1600
1601cpl_error_code gravi_acqcam_get_pupil_offset_v2(cpl_imagelist ** pupilImage_shiftandadd,
1602 cpl_array * bad_frames_short,
1603 cpl_bivector ** diode_pos_offset, cpl_propertylist * o_header,
1604 cpl_size nrow_on)
1605{
1607
1608 cpl_ensure_code(pupilImage_shiftandadd, CPL_ERROR_NULL_INPUT);
1609 cpl_ensure_code(bad_frames_short, CPL_ERROR_NULL_INPUT);
1610 cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1611 cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
1612 cpl_ensure_code(nrow_on, CPL_ERROR_NULL_INPUT);
1613
1614 int nv =0;
1615
1616 if (cpl_array_get_size(bad_frames_short)!= nrow_on)
1617 cpl_msg_error (cpl_func, "Problem with number of blinking frames");
1618
1619 for (int tel = 0 ; tel < GRAVI_SPOT_NTEL; tel++)
1620 {
1621 double previous_xsc=0.0;
1622 double previous_ysc=0.0;
1623 double flux_sum=0.0;
1624 int n_sum=0;
1625
1626 for (cpl_size n = 0 ; n < nrow_on; n++)
1627 {
1628 cpl_image * image = cpl_imagelist_get (pupilImage_shiftandadd[tel], n);
1629 CPLCHECK_MSG("Failure at getting saa image");
1630
1631 cpl_size px, py;
1632 cpl_image_get_maxpos (image, &px, &py);
1633 CPLCHECK_MSG("Failure to read maximum position");
1634
1635 double flux_max= cpl_image_get(image,px,py,&nv);
1636 double std = cpl_image_get_stdev (image);
1637
1638
1639 double xsc = (double) px;
1640 double ysc = (double) py;
1641 double exsc=3., eysc=3.;
1642 cpl_size size = 15;
1643
1644 gravi_acq_fit_gaussian (image, &xsc, &ysc, &exsc, &eysc, size);
1645 CPLCHECK_MSG("Failure at fitting pupil spot");
1646
1647 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n]);
1648 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n]);
1649 CPLCHECK_MSG("Failure at reading bivector offsets");
1650
1651 double x_offset_final= cpl_vector_get(x_pos_offset, tel) + xsc-GRAVI_SPOT_NSEARCH-1;
1652 double y_offset_final= cpl_vector_get(y_pos_offset, tel) + ysc-GRAVI_SPOT_NSEARCH-1;
1653
1654 cpl_vector_set(x_pos_offset, tel, x_offset_final);
1655 cpl_vector_set(y_pos_offset, tel, y_offset_final);
1656
1657 int previous;
1658 if ((flux_max/std < 4)||(exsc<0))
1659 {
1660 previous = cpl_array_get_int (bad_frames_short,n,&nv);
1661 cpl_array_set_int (bad_frames_short,n,previous|(1<<tel));
1662 } else {
1663 if (n>0)
1664 {
1665 double distance=(xsc-previous_xsc)*(xsc-previous_xsc)+(ysc-previous_ysc)*(ysc-previous_ysc);
1666 if (distance > 10*10)
1667 {
1668 previous = cpl_array_get_int (bad_frames_short,n-1,&nv);
1669 cpl_array_set_int (bad_frames_short,n-1,previous|(1<<(tel+4)));
1670 previous = cpl_array_get_int (bad_frames_short,n,&nv);
1671 cpl_array_set_int (bad_frames_short,n,previous|(1<<(tel+4)));
1672 }
1673 }
1674 flux_sum+=flux_max;
1675 n_sum+=1;
1676 }
1677 previous_xsc=xsc;
1678 previous_ysc=ysc;
1679 CPLCHECK_MSG("Failure at testing bad pupil measurments");
1680 }
1681
1682 char qc_name[100];
1683 double flux_mean=flux_sum/(1e-7+n_sum);
1684 /* add QC parameters with pupil beacon flux */
1685 sprintf(qc_name, "ESO QC ACQ PUP%i SPOT_FLUX", tel + 1);
1686 cpl_msg_info(cpl_func, "%s = %f", qc_name,flux_mean);
1687 cpl_propertylist_update_double(o_header, qc_name, flux_mean);
1688 cpl_propertylist_set_comment(o_header, qc_name,"[ADU] mean spot flux");
1689
1690 CPLCHECK_MSG("Failing to store QC value");
1691 }
1692
1693 CPLCHECK_MSG("Failure at fitting gaussian on diode position");
1694
1696 return CPL_ERROR_NONE;
1697}
1698
1699
1700/*----------------------------------------------------------------------------*/
1722/*----------------------------------------------------------------------------*/
1723
1724
1725
1726cpl_error_code gravi_acqcam_set_pupil_table_v2(cpl_table * acqcam_table, cpl_propertylist * header, cpl_propertylist * o_header, cpl_array * good_frames,
1727 cpl_vector* scale_vector, cpl_array * bad_frames_short, cpl_bivector ** diode_pos_offset, cpl_vector * focus_value, gravi_data *static_param_data)
1728{
1729
1731int nv =0;
1732
1733cpl_ensure_code(acqcam_table, CPL_ERROR_NULL_INPUT);
1734cpl_ensure_code(header, CPL_ERROR_NULL_INPUT);
1735cpl_ensure_code(o_header, CPL_ERROR_NULL_INPUT);
1736cpl_ensure_code(good_frames, CPL_ERROR_NULL_INPUT);
1737cpl_size nrow = cpl_array_get_size(good_frames);
1738
1739/* Pupil positions array */
1740gravi_table_new_column(acqcam_table, "PUPIL_NSPOT", NULL, CPL_TYPE_INT);
1741gravi_table_new_column(acqcam_table, "PUPIL_X", "pix", CPL_TYPE_DOUBLE);
1742gravi_table_new_column(acqcam_table, "PUPIL_Y", "pix", CPL_TYPE_DOUBLE);
1743gravi_table_new_column(acqcam_table, "PUPIL_Z", "pix", CPL_TYPE_DOUBLE);
1744gravi_table_new_column(acqcam_table, "PUPIL_R", "deg", CPL_TYPE_DOUBLE);
1745gravi_table_new_column(acqcam_table, "PUPIL_U", "m", CPL_TYPE_DOUBLE);
1746gravi_table_new_column(acqcam_table, "PUPIL_V", "m", CPL_TYPE_DOUBLE);
1747gravi_table_new_column(acqcam_table, "PUPIL_W", "m", CPL_TYPE_DOUBLE);
1748gravi_table_new_column(acqcam_table, "OPD_PUPIL", "m", CPL_TYPE_DOUBLE);
1749
1750double sobj_x = gravi_pfits_get_sobj_x(header);
1751double sobj_y = gravi_pfits_get_sobj_y(header);
1752CPLCHECK_MSG("Cannot determine SOBJ X and Y");
1753
1754if (cpl_array_get_mean(good_frames)*nrow<.5)
1755{
1756 /* no good frames */
1757 /* filling the table with zeros */
1758 cpl_msg_info(cpl_func,"No good pupil images ==> no good pupil position available");
1759
1760 for (int tel = 0; tel < GRAVI_SPOT_NTEL; tel++)
1761 for (cpl_size row = 0; row < nrow; row++)
1762 {
1763 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1764 CPLCHECK_MSG("Cannot put 0 in NSPOT data in the ACQ PUPIL");
1765 }
1766
1768 return CPL_ERROR_NONE;
1769}
1770
1771 /* looks like the table will not be empty. Checking additional parameters*/
1772cpl_ensure_code(scale_vector, CPL_ERROR_NULL_INPUT);
1773cpl_ensure_code(bad_frames_short, CPL_ERROR_NULL_INPUT);
1774cpl_ensure_code(diode_pos_offset, CPL_ERROR_NULL_INPUT);
1775cpl_ensure_code(focus_value, CPL_ERROR_NULL_INPUT);
1776cpl_ensure_code(static_param_data, CPL_ERROR_NULL_INPUT);
1777
1778cpl_size nrow_on = cpl_array_get_size(bad_frames_short);
1779
1780if (fabs(cpl_array_get_mean(good_frames)*nrow -nrow_on) > 1e-2)
1781 cpl_msg_error (cpl_func, "Ratio of blinking frames different %f from %lli",cpl_array_get_mean(good_frames)*nrow,nrow_on);
1782
1783
1784
1785for (int tel = 0; tel < GRAVI_SPOT_NTEL; tel++)
1786{
1787 int n_on=0;
1788
1789 /* Get the conversion angle xy to uv in [rad] */
1790 double northangle = gravi_pfits_get_northangle_acqcam(header, tel);
1791 double cfangle = cos(northangle * CPL_MATH_RAD_DEG);
1792 double sfangle = sin(northangle * CPL_MATH_RAD_DEG);
1793 double scale = cpl_vector_get(scale_vector,tel);
1794 CPLCHECK_MSG("Cannot determine field angle");
1795
1796 /* values which are valids for all frames */
1797 double r_shift = 0.0;
1798 double z_shift = cpl_vector_get(focus_value,tel)/GRAVI_SPOT_SWINDOW;
1799 double w_shift = gravi_acqcam_z2meter(z_shift, static_param_data);
1800
1801 /*values that does change */
1802 double x_shift = 0.0;
1803 double y_shift = 0.0;
1804 double u_shift = 0.0;
1805 double v_shift = 0.0;
1806 double opd_pupil = 0.0;
1807
1808 /* buffers to average values */
1809 double x_shift_sum = 0.0;
1810 double y_shift_sum = 0.0;
1811 double u_shift_sum = 0.0;
1812 double v_shift_sum = 0.0;
1813 int n_sum = 0;
1814 int n_bad_snr = 0;
1815 int n_bad_distance = 0;
1816
1817 for (cpl_size row = 0; row < nrow; row++)
1818 {
1819 /*cpl_msg_info (cpl_func, "Row number %lli n %lli",row,n_on);*/
1820
1821 if (cpl_array_get_int(good_frames,row,&nv) != 1)
1822 {
1823 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1824 CPLCHECK_MSG("Cannot put 0 in NSPOT data in the ACQ PUPIL");
1825 }
1826 else
1827 {
1828
1829 int is_it_bad=cpl_array_get_int (bad_frames_short,n_on,&nv);
1830
1831 cpl_vector * x_pos_offset = cpl_bivector_get_x (diode_pos_offset[n_on]);
1832 cpl_vector * y_pos_offset = cpl_bivector_get_y (diode_pos_offset[n_on]);
1833 CPLCHECK_MSG("Cannot put data from the offset bivector");
1834
1835 x_shift = cpl_vector_get(x_pos_offset,tel);
1836 y_shift = cpl_vector_get(y_pos_offset,tel);
1837
1838 /* In UV [m] */
1839 u_shift = (cfangle * x_shift - sfangle * y_shift)
1840 / scale;
1841 v_shift = (sfangle * x_shift + cfangle * y_shift)
1842 / scale;
1843 opd_pupil = -(u_shift * sobj_x + v_shift * sobj_y)
1845
1846 CPLCHECK_MSG("Cannot prepare data to be put in the ACQ PUPIL table");
1847
1848
1849 if (CHECK_BIT(is_it_bad,tel)||CHECK_BIT(is_it_bad,tel+4))
1850 {
1851 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel, 0);
1852
1853 /* increasing counters */
1854 if (CHECK_BIT(is_it_bad,tel )) n_bad_snr += 1;
1855 if (CHECK_BIT(is_it_bad,tel+4)) n_bad_distance += 1;
1856 u_shift = 0.0;
1857 v_shift = 0.0;
1858 opd_pupil = 0.0;
1859 cpl_msg_info(cpl_func,"Pupil image number%4lli is un-usable for tel %i",row,tel);
1860 }
1861 else
1862 {
1863 cpl_table_set(acqcam_table, "PUPIL_NSPOT", row * GRAVI_SPOT_NTEL + tel,
1864 16);
1865
1866 /* averaging positions */
1867 x_shift_sum += x_shift;
1868 y_shift_sum += y_shift;
1869 u_shift_sum += u_shift;
1870 v_shift_sum += v_shift;
1871
1872 /* increasing counter */
1873 n_sum += 1;
1874 }
1875
1876 cpl_table_set(acqcam_table, "PUPIL_X", row * GRAVI_SPOT_NTEL + tel,
1877 x_shift);
1878 cpl_table_set(acqcam_table, "PUPIL_Y", row * GRAVI_SPOT_NTEL + tel,
1879 y_shift);
1880 cpl_table_set(acqcam_table, "PUPIL_Z", row * GRAVI_SPOT_NTEL + tel,
1881 z_shift);
1882 cpl_table_set(acqcam_table, "PUPIL_R", row * GRAVI_SPOT_NTEL + tel,
1883 r_shift);
1884 cpl_table_set(acqcam_table, "PUPIL_U", row * GRAVI_SPOT_NTEL + tel,
1885 u_shift);
1886 cpl_table_set(acqcam_table, "PUPIL_V", row * GRAVI_SPOT_NTEL + tel,
1887 v_shift);
1888 cpl_table_set(acqcam_table, "PUPIL_W", row * GRAVI_SPOT_NTEL + tel,
1889 w_shift);
1890 cpl_table_set(acqcam_table, "OPD_PUPIL", row * GRAVI_SPOT_NTEL + tel,
1891 opd_pupil);
1892 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1893
1894 n_on++;
1895 }
1896 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1897 }
1898
1899 /* Add QC parameters */
1900 char qc_name[100];
1901
1902 sprintf(qc_name, "ESO QC ACQ FIELD%i NORTH_ANGLE", tel + 1);
1903 cpl_msg_info(cpl_func, "%s = %f", qc_name, northangle);
1904 cpl_propertylist_update_double(o_header, qc_name, northangle);
1905 cpl_propertylist_set_comment(o_header, qc_name,
1906 "[deg] y->x, predicted North direction on ACQ");
1907
1908 sprintf(qc_name, "ESO QC ACQ PUP%i FRAMES", tel + 1);
1909 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_sum);
1910 cpl_propertylist_update_double(o_header, qc_name, n_sum);
1911 cpl_propertylist_set_comment(o_header, qc_name,
1912 "Good ACQ pupil frames");
1913
1914 sprintf(qc_name, "ESO QC ACQ PUP%i BADSNR", tel + 1);
1915 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_bad_snr);
1916 cpl_propertylist_update_double(o_header, qc_name, n_bad_snr);
1917 cpl_propertylist_set_comment(o_header, qc_name,
1918 "Frames with low SNR");
1919
1920 sprintf(qc_name, "ESO QC ACQ PUP%i JUMP", tel + 1);
1921 cpl_msg_info(cpl_func, "%s = %i", qc_name, n_bad_distance);
1922 cpl_propertylist_update_double(o_header, qc_name, n_bad_distance);
1923 cpl_propertylist_set_comment(o_header, qc_name,
1924 "Frames with position jump");
1925
1926 sprintf(qc_name, "ESO QC ACQ PUP%i SCALE", tel + 1);
1927 cpl_msg_info(cpl_func, "%s = %f", qc_name, scale);
1928 cpl_propertylist_update_double(o_header, qc_name, scale);
1929 cpl_propertylist_set_comment(o_header, qc_name,
1930 "[pix/m] diode scale on ACQ");
1931
1932 sprintf(qc_name, "ESO QC ACQ PUP%i FOCUS", tel + 1);
1933 cpl_msg_info(cpl_func, "%s = %f", qc_name, z_shift);
1934 cpl_propertylist_update_int(o_header, qc_name, z_shift);
1935 cpl_propertylist_set_comment(o_header, qc_name,
1936 "[pix] defocus of pupil plane");
1937
1938 /* to avoid division by zero */
1939 if (n_sum ==0) n_sum+=1;
1940
1941 sprintf(qc_name, "ESO QC ACQ PUP%i XPOS", tel + 1);
1942 cpl_msg_info(cpl_func, "%s = %f", qc_name, x_shift_sum/n_sum);
1943 cpl_propertylist_update_double(o_header, qc_name, x_shift_sum/n_sum);
1944 cpl_propertylist_set_comment(o_header, qc_name,
1945 "[pix] pupil x-shift in ACQ");
1946
1947 sprintf(qc_name, "ESO QC ACQ PUP%i YPOS", tel + 1);
1948 cpl_msg_info(cpl_func, "%s = %f", qc_name, y_shift_sum/n_sum);
1949 cpl_propertylist_update_double(o_header, qc_name, y_shift_sum/n_sum);
1950 cpl_propertylist_set_comment(o_header, qc_name,
1951 "[pix] pupil y-shift in ACQ");
1952
1953 sprintf(qc_name, "ESO QC ACQ PUP%i UPOS", tel + 1);
1954 cpl_msg_info(cpl_func, "%s = %f", qc_name, u_shift_sum/n_sum);
1955 cpl_propertylist_update_double(o_header, qc_name, u_shift_sum/n_sum);
1956 cpl_propertylist_set_comment(o_header, qc_name,
1957 "[m] pupil u-shift in ACQ");
1958
1959 sprintf(qc_name, "ESO QC ACQ PUP%i VPOS", tel + 1);
1960 cpl_msg_info(cpl_func, "%s = %f", qc_name, v_shift_sum/n_sum);
1961 cpl_propertylist_update_double(o_header, qc_name, v_shift_sum/n_sum);
1962 cpl_propertylist_set_comment(o_header, qc_name,
1963 "[m] pupil v-shift in ACQ");
1964
1965
1966 CPLCHECK_MSG("Cannot put data in the ACQ PUPIL table");
1967
1968
1969
1970}
1971
1973return CPL_ERROR_NONE;
1974}
1975
1976
1977
1978/*----------------------------------------------------------------------------*/
1995/*----------------------------------------------------------------------------*/
1996
1997
1998
1999cpl_image * gravi_image_extract(cpl_image * image_in, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
2000 {
2001
2002 cpl_size nx = cpl_image_get_size_x (image_in);
2003 cpl_size ny = cpl_image_get_size_y (image_in);
2004 cpl_size llx_new,lly_new,urx_new,ury_new;
2005 cpl_size xpos, ypos;
2006
2007 cpl_size nx_new=urx+1-llx;
2008 cpl_size ny_new=ury+1-lly;
2009 cpl_image * image_out= cpl_image_new (nx_new, ny_new, cpl_image_get_type(image_in));
2010 cpl_image_fill_window (image_out, 1, 1, nx_new, ny_new, 0.0);
2011
2012
2013 if ((llx >= nx)||(lly >= ny)||(urx <= 0)||(ury <= 0))
2014 {
2015 cpl_msg_warning (cpl_func, "Cutting at x=(%lli,%lli) y=(%lli,%lli) is outside the window bondaries", llx, urx, lly, ury);
2016 }
2017 else
2018 {
2019 if (llx < 1)
2020 {
2021 llx_new=1;
2022 xpos=1-llx;
2023 } else
2024 {
2025 llx_new=llx;
2026 xpos=1;
2027 }
2028 if (lly < 1)
2029 {
2030 lly_new=1;
2031 ypos=1-lly;
2032 } else
2033 {
2034 lly_new=lly;
2035 ypos=1;
2036 }
2037
2038 if (urx > nx)
2039 urx_new=nx;
2040 else
2041 urx_new=urx;
2042 if (ury > ny)
2043 ury_new=ny;
2044 else
2045 ury_new=ury;
2046
2047 cpl_image * image_in_cut = cpl_image_extract(image_in,llx_new,lly_new,urx_new,ury_new);
2048
2049 cpl_image_copy (image_out, image_in_cut, xpos, ypos);
2050
2051 cpl_image_delete(image_in_cut);
2052 }
2053
2054 return image_out;
2055 }
2056
2057/*----------------------------------------------------------------------------*/
2062/*----------------------------------------------------------------------------*/
2063
2064
2066{
2067 double defocus;
2068 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);
2069
2070 return defocus;
2071}
2072
2073/*----------------------------------------------------------------------------*/
2090/*----------------------------------------------------------------------------*/
2091
2092cpl_error_code gravi_acqcam_field (cpl_image * mean_img,
2093 cpl_imagelist * acqcam_imglist,
2094 cpl_propertylist * header,
2095 cpl_table * acqcam_table,
2096 cpl_propertylist * o_header,
2097 gravi_data *static_param_data)
2098{
2100 cpl_ensure_code (mean_img, CPL_ERROR_NULL_INPUT);
2101 cpl_ensure_code (acqcam_imglist, CPL_ERROR_NULL_INPUT);
2102 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
2103 cpl_ensure_code (acqcam_table, CPL_ERROR_NULL_INPUT);
2104 cpl_ensure_code (o_header, CPL_ERROR_NULL_INPUT);
2105
2106 char qc_name[100];
2107 int ntel = 4;
2108
2109 int sts_mode = 0;
2110
2111 /* check the feed mode */
2112 if ( !strcmp(gravi_pfits_get_feed (header), "DUAL_STS") ) {
2113 sts_mode = DUAL_STS;
2114 }
2115 else{
2116 sts_mode = SINGLE_STS;
2117 }
2118
2119 /* Number of row */
2120 cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
2121
2122 /* MODE_ONAXIS or MODE_OFFAIS. We query once
2123 at the begining of the function */
2124 int axis_mode = gravi_pfits_get_axis (header);
2125
2126 /* Create columns */
2127 gravi_table_new_column (acqcam_table, "FIELD_SC_X", "pix", CPL_TYPE_DOUBLE);
2128 gravi_table_new_column (acqcam_table, "FIELD_SC_Y", "pix", CPL_TYPE_DOUBLE);
2129 gravi_table_new_column (acqcam_table, "FIELD_SC_XERR", "pix", CPL_TYPE_DOUBLE);
2130 gravi_table_new_column (acqcam_table, "FIELD_SC_YERR", "pix", CPL_TYPE_DOUBLE);
2131 gravi_table_new_column (acqcam_table, "FIELD_FT_X", "pix", CPL_TYPE_DOUBLE);
2132 gravi_table_new_column (acqcam_table, "FIELD_FT_Y", "pix", CPL_TYPE_DOUBLE);
2133 gravi_table_new_column (acqcam_table, "FIELD_FT_XERR", "pix", CPL_TYPE_DOUBLE);
2134 gravi_table_new_column (acqcam_table, "FIELD_FT_YERR", "pix", CPL_TYPE_DOUBLE);
2135 gravi_table_new_column (acqcam_table, "FIELD_SCALE", "mas/pix", CPL_TYPE_DOUBLE);
2136 gravi_table_new_column (acqcam_table, "FIELD_SCALEERR", "mas/pix", CPL_TYPE_DOUBLE);
2137 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DX", "pix", CPL_TYPE_DOUBLE);
2138 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DY", "pix", CPL_TYPE_DOUBLE);
2139 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DXERR", "pix", CPL_TYPE_DOUBLE);
2140 gravi_table_new_column (acqcam_table, "FIELD_FIBER_DYERR", "pix", CPL_TYPE_DOUBLE);
2141 gravi_table_new_column (acqcam_table, "FIELD_STREHL", "ratio", CPL_TYPE_DOUBLE);
2142
2143 /* ----------------- START EKW 26/11/2018 read constant parameter from calibration file */
2144 /* Position of roof center on full frame */
2145 //double roof_x[] = {274.4, 787.1, 1236.1, 1673.4};
2146 //double roof_y[] = {242.3, 247.7, 225.8, 235.6};
2147 double *roof_x = NULL;
2148 double *roof_y = NULL;
2149
2150 /* Position of single-field spot on full frame */
2151 //double spot_x[] = {289. , 798.2, 1245.5, 1696.};
2152 //double spot_y[] = {186.5, 187.5, 172.5, 178.};
2153 double *spot_x = NULL;
2154 double *spot_y = NULL;
2155
2156 /* Default position angle of roof */
2157 //double roof_pos[] = {38.49, 38.54, 38.76, 39.80};
2158 double *roof_pos = NULL;
2159
2160 cpl_table * roof_pos_table = gravi_data_get_table (static_param_data, "ROOFPOS");
2161 CPLCHECK_MSG ("STATIC_PARAM not available in the SOF. It is mandatory for acqcam reduction.");
2162
2163 if ( cpl_table_has_column(roof_pos_table , "roof_x") ) {
2164 roof_x= cpl_table_get_data_double (roof_pos_table, "roof_x");
2165 cpl_msg_info(cpl_func,"roof_x [0] : %e", roof_x[0] );
2166 cpl_msg_info(cpl_func,"roof_x [1] : %e", roof_x[1] );
2167 cpl_msg_info(cpl_func,"roof_x [2] : %e", roof_x[2] );
2168 cpl_msg_info(cpl_func,"roof_x [3] : %e", roof_x[3] );
2169 }
2170 else {
2171 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_x ");
2172 }
2173
2174 if ( cpl_table_has_column(roof_pos_table , "roof_y") ) {
2175 roof_y= cpl_table_get_data_double (roof_pos_table, "roof_y");
2176 cpl_msg_info(cpl_func,"roof_y [0] : %e", roof_y[0] );
2177 cpl_msg_info(cpl_func,"roof_y [1] : %e", roof_y[1] );
2178 cpl_msg_info(cpl_func,"roof_y [2] : %e", roof_y[2] );
2179 cpl_msg_info(cpl_func,"roof_y [3] : %e", roof_y[3] );
2180 }
2181 else {
2182 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_y ");
2183 }
2184
2185 if ( cpl_table_has_column(roof_pos_table , "spot_x") ) {
2186 spot_x= cpl_table_get_data_double (roof_pos_table, "spot_x");
2187 cpl_msg_info(cpl_func,"spot_x [0] : %e", spot_x[0] );
2188 cpl_msg_info(cpl_func,"spot_x [1] : %e", spot_x[1] );
2189 cpl_msg_info(cpl_func,"spot_x [2] : %e", spot_x[2] );
2190 cpl_msg_info(cpl_func,"spot_x [3] : %e", spot_x[3] );
2191 }
2192 else {
2193 cpl_msg_warning(cpl_func,"Cannot get the default values for spot_x ");
2194 }
2195
2196 if ( cpl_table_has_column(roof_pos_table , "spot_y") ) {
2197 spot_y= cpl_table_get_data_double (roof_pos_table, "spot_y");
2198 cpl_msg_info(cpl_func,"spot_y [0] : %e", spot_y[0] );
2199 cpl_msg_info(cpl_func,"spot_y [1] : %e", spot_y[1] );
2200 cpl_msg_info(cpl_func,"spot_y [2] : %e", spot_y[2] );
2201 cpl_msg_info(cpl_func,"spot_y [3] : %e", spot_y[3] );
2202 }
2203 else {
2204 cpl_msg_warning(cpl_func,"Cannot get the default values for spot_y ");
2205 }
2206
2207 if ( cpl_table_has_column(roof_pos_table , "roof_pos") ) {
2208 roof_pos= cpl_table_get_data_double (roof_pos_table, "roof_pos");
2209 cpl_msg_info(cpl_func,"roof_pos [0] : %e", roof_pos[0] );
2210 cpl_msg_info(cpl_func,"roof_pos [1] : %e", roof_pos[1] );
2211 cpl_msg_info(cpl_func,"roof_pos [2] : %e", roof_pos[2] );
2212 cpl_msg_info(cpl_func,"roof_pos [3] : %e", roof_pos[3] );
2213 }
2214 else {
2215 cpl_msg_warning(cpl_func,"Cannot get the default values for roof_pos ");
2216 }
2217 /* ------------------ END EKW 26/11/2018 read constant parameter from calibration file */
2218
2219
2220
2221 /* Static acqcam field calibration error */
2222 /* calibration from Frank, March 2017, using Marcel
2223 double calib_dx[] = {-0.3483, -1.0251, -0.5432, -0.2024} ;
2224 double calib_dy[] = { 0.3089, -0.5221, -0.2686, -0.3843} ; */
2225 /* calibration from Julien, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-10
2226 double calib_dx[] = {-0.61, -1.69, -0.97, 0.00} ;
2227 double calib_dy[] = { 0.53, -1.07, -0.49, -0.47} ; */
2228 /* calibration from Oli, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-19/
2229 double calib_dx[] = {-0.63, -1.67, -0.97, 0.02} ;
2230 double calib_dy[] = { 0.55, -1.07, -0.49, -0.45} ; */
2231 /* FE: removing obsolete calib_dx, because now handled by ICS, or by
2232 correcting fitsheader retrospectively */
2233 double calib_dx[] = {0,0,0,0} ;
2234 double calib_dy[] = {0,0,0,0} ;
2235
2236 /* If sub-windowing, we read the sub-window size */
2237 cpl_size nsx = 512;
2238 /* EKW 10/01/2019 cpl_size nsy = 512; */
2239 if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
2240 nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
2241 /* EKW 10/01/2019 nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY"); */
2242 }
2243
2244 /* Compute separation */
2245 double sobj_x = gravi_pfits_get_sobj_x (header);
2246 double sobj_y = gravi_pfits_get_sobj_y (header);
2247
2248 /* Accumulated mapping offsets */
2249 double sobj_offx=0., sobj_offy=0, sobj_drho=0., sobj_dth=0.;
2250 if (cpl_propertylist_has(header, "ESO INS SOBJ OFFX")) {
2251 /* mapping / mosaicing blind offsets from otiginal acquisition
2252 * to current position */
2253 sobj_offx = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFX");
2254 sobj_offy = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFY");
2255
2256 /* distance from acquired SC position to current SC position in mas */
2257 sobj_drho = sqrt(sobj_offx*sobj_offx+sobj_offy*sobj_offy);
2258
2259 /* position angle on sky from acquired SC position to */
2260 /* current SC position, neglecting anamorphism variations */
2261 sobj_dth = atan2(sobj_offx, sobj_offy)*180./M_PI;
2262 }
2263 CPLCHECK_MSG ("Cannot get separation");
2264
2265 /* Recover position of originally acquired star, before any blind offset */
2266 double dx_in = sobj_x - sobj_offx;
2267 double dy_in = sobj_y - sobj_offy;
2268 double rho_in = sqrt(dx_in*dx_in + dy_in*dy_in);
2269
2270 /* Force zero separation for SINGLE mode */
2271 if ( gravi_pfits_get_mode (header) == MODE_SINGLE) rho_in = 0.;
2272 CPLCHECK_MSG ("Error reading header information");
2273
2274 /* Loop on tel */
2275 for (int tel = 0; tel < ntel; tel++) {
2276 char name[90];
2277 double rp=roof_pos[tel]; // default value of roof position angle
2278 cpl_size sx=tel*512+1, sy=1;
2279 cpl_msg_info (cpl_func, "Compute field position for beam %i", tel+1);
2280
2281 /* Read roof position angle */
2282 sprintf (name, "ESO INS DROTOFF%d", tel + 1);
2283 if ( cpl_propertylist_has (header, name) ) {
2284 rp = cpl_propertylist_get_double(header, name);
2285 CPLCHECK ("Cannot get rotation");
2286 }
2287 else {
2288 cpl_msg_info (cpl_func, "%s not in header, use %f", name, rp);
2289 }
2290
2291 /* Approx. position angle of the binary, left from top */
2292 double approx_PA = 270.-rp;
2293
2294 /* Get the telescope name and ID */
2295 const char * telname = gravi_conf_get_telname (tel, header);
2296
2297 /* Check telescope name */
2298 if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
2299 cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
2300
2301 /* Hardcoded approx. plate-scale in mas/pix. */
2302 double scale = 0.0;
2303 char telid[128];
2304 if (telname[0] == 'U') {
2305 /* scale = 18.; */
2306 scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO PLATE SCALE UT");
2307 cpl_msg_info (cpl_func,"PLATE SCALE UT is : %e",scale);
2308 } else if (telname[0] == 'A') {
2309 sprintf(telid, "ESO PLATE SCALE AT%d", tel+1);
2310 scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), telid);
2311 cpl_msg_info (cpl_func,"PLATE SCALE AT%d is : %e",tel+1, scale);
2312 /*
2313 if (tel == 0) {
2314 scale = 76.8;
2315 } else if (tel == 1) {
2316 scale = 78.0;
2317 } else if (tel == 2) {
2318 scale = 77.0;
2319 } else if (tel == 3) {
2320 scale = 84.6;
2321 }
2322 */
2323 }
2324
2325 /* Position of the fibres */
2326 double fiber_xft=0.;
2327 double fiber_yft=0.;
2328 double fiber_xsc=0.;
2329 double fiber_ysc=0.;
2330 sprintf (name, "ESO ACQ FIBER FT%dX", tel + 1);
2331 if (cpl_propertylist_has(header, name)) {
2332 fiber_xft=cpl_propertylist_get_double(header, name);
2333 sprintf (name, "ESO ACQ FIBER FT%dY", tel + 1);
2334 fiber_yft=cpl_propertylist_get_double(header, name);
2335 sprintf (name, "ESO ACQ FIBER SC%dX", tel + 1);
2336 fiber_xsc=cpl_propertylist_get_double(header, name);
2337 sprintf (name, "ESO ACQ FIBER SC%dY", tel + 1);
2338 fiber_ysc=cpl_propertylist_get_double(header, name);
2339 }
2340 double fiber_ft_sc_x=fiber_xsc-fiber_xft;
2341 double fiber_ft_sc_y=fiber_ysc-fiber_yft;
2342
2343 /* Get the North position angle on the camera */
2344 double northangle = gravi_pfits_get_northangle_acqcam (header, tel);
2345 CPLCHECK ("Cannot get rotation");
2346
2347 /* Mapping/mosaicing offset on acq cam axes, in mas, */
2348 /* neglecting amnamorphism variations */
2349 double sobj_offx_cam = sobj_drho * sin((northangle+sobj_dth)/180.*M_PI);
2350 double sobj_offy_cam = sobj_drho * cos((northangle+sobj_dth)/180.*M_PI);
2351
2352 /* If sub-windowing, we read the sub-window start for field */
2353 if ( nsx != 512 ) {
2354 sprintf (name, "ESO DET1 FRAM%d STRX", tel + 1);
2355 sx = cpl_propertylist_get_int (header, name);
2356
2357 sprintf (name, "ESO DET1 FRAM%d STRY", tel + 1);
2358 sy = cpl_propertylist_get_int (header, name);
2359
2360 CPLCHECK_MSG ("Cannot get sub-windowing parameters");
2361 }
2362
2363 cpl_msg_debug (cpl_func,"sub-window field %i sx= %lld sy = %lld", tel, sx, sy);
2364
2365 /*--------------------------------------------------*/
2366 /* Guess position of SC and FT targes in mean image */
2367 /*--------------------------------------------------*/
2368
2369 double xFT, yFT, xSC, ySC;
2370
2371 if (axis_mode == MODE_ONAXIS) {
2372 /* TODO: close dual-field */
2373 /* Single-field case */
2374 /* Simply shift the best spot from full frame to cut-out */
2375 xFT = spot_x[tel] - sx + nsx*tel + 1;
2376 yFT = spot_y[tel] - sy + 1;
2377 xSC = xFT;
2378 ySC = yFT;
2379 } else if (sts_mode == DUAL_STS) {
2380 xFT = fiber_xft - sx + nsx*tel + 1;
2381 yFT = fiber_yft - sy + 1;
2382 xSC = fiber_xsc - sx + nsx*tel + 1;
2383 ySC = fiber_ysc - sy + 1;
2384 } else {
2385 /* Pixel position of roof center on cut-out frame */
2386 /* Shift from full frame to cut-out */
2387 double cutout_roof_x = roof_x[tel] - sx + nsx*tel + 1;
2388 double cutout_roof_y = roof_y[tel] - sy + 1;
2389
2390 cpl_msg_info (cpl_func, " ROOF_X = %f", cutout_roof_x);
2391 cpl_msg_info (cpl_func, " ROOF_Y = %f", cutout_roof_y);
2392
2393 /* Approx pixel offset from SC to FT, divided by 2 */
2394 double approx_dx=0.5*rho_in*sin(approx_PA*M_PI/180.)/scale;
2395 double approx_dy=0.5*rho_in*cos(approx_PA*M_PI/180.)/scale;
2396
2397 /* Expected position of the two stars */
2398 xFT = cutout_roof_x - approx_dx ;
2399 yFT = cutout_roof_y - approx_dy ;
2400 xSC = cutout_roof_x + approx_dx ;
2401 ySC = cutout_roof_y + approx_dy ;
2402 }
2403
2404 cpl_msg_info (cpl_func, "guess SC_X = %f", xSC);
2405 cpl_msg_info (cpl_func, "guess SC_Y = %f", ySC);
2406 cpl_msg_info (cpl_func, "guess FT_X = %f", xFT);
2407 cpl_msg_info (cpl_func, "guess FT_Y = %f", yFT);
2408
2409 /*------------------------------------------*/
2410 /* Measure SC target position in mean image */
2411 /*------------------------------------------*/
2412
2413 /* Box size */
2414 /* Optimal size has been determined empirically. */
2415 /* Too small value (15) will sometimes miss even a bright target. */
2416 /* Too large value (25) will often pick a wrong star (or backround noise) */
2417 cpl_size size = 20;
2418
2419 double xSCguess=xSC, ySCguess=ySC, xFTguess=xFT, yFTguess=yFT;
2420 double ex, ey;
2421 double qc_val=0.;
2422
2423 gravi_acq_fit_gaussian (mean_img, &xSC, &ySC, &ex, &ey, size);
2424 CPLCHECK_MSG("Error fitting SC");
2425
2426 /*--------------------------------------------------------*/
2427 /* Add QC parameters for SC target position in mean image */
2428 /*--------------------------------------------------------*/
2429
2430 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X", tel+1);
2431 if (xSC==0.) {
2432 /* Fitting failed: put QC to 0., reset xSC to gues value */
2433 qc_val = 0.;
2434 xSC = xSCguess;
2435 } else {
2436 /* Fiting succeeded: shift into full frame */
2437 qc_val = xSC;
2438 }
2439 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2440 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2441 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2442
2443 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y", tel+1);
2444 if (ySC==0.) {
2445 /* Fitting failed: put QC to 0., reset ySC to gues value */
2446 qc_val = 0.;
2447 ySC = ySCguess;
2448 } else {
2449 /* Fiting succeeded: shift into full frame */
2450 qc_val = ySC;
2451 }
2452 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2453 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2454 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2455
2456 /*------------------------------------------*/
2457 /* Measure FT target position in mean image */
2458 /*------------------------------------------*/
2459
2460 if (axis_mode != MODE_ONAXIS) {
2461 /* Detec FT */
2462 gravi_acq_fit_gaussian (mean_img, &xFT, &yFT, &ex, &ey, size);
2463 CPLCHECK_MSG("Error fitting FT");
2464 } else {
2465 xFT=xSC;
2466 yFT=ySC;
2467 }
2468
2469 /*--------------------------------------------------------*/
2470 /* Add QC parameters for FT target position in mean image */
2471 /*--------------------------------------------------------*/
2472
2473 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X", tel+1);
2474 if (xFT==0.) {
2475 /* Fitting failed: put QC to 0., reset xFT to gues value */
2476 qc_val = 0.;
2477 xFT = xFTguess;
2478 } else {
2479 /* Fiting succeeded: shift into full frame */
2480 qc_val = xFT;
2481 }
2482 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2483 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2484 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2485
2486 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y", tel+1);
2487 if (yFT==0.) {
2488 /* Fitting failed: put QC to 0., reset yFT to gues value */
2489 qc_val = 0.;
2490 yFT = yFTguess;
2491 } else {
2492 /* Fiting succeeded: shift into full frame */
2493 qc_val = yFT;
2494 }
2495 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2496 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2497 cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
2498
2499 /*---------------------------------------*/
2500 /* If in dual field, measure plate scale */
2501 /*---------------------------------------*/
2502
2503 if (axis_mode != MODE_ONAXIS) {
2504 sprintf (qc_name, "ESO QC ACQ FIELD%i SCALE", tel+1);
2505 double sep = sqrt((ySC-yFT)*(ySC-yFT)+(xSC-xFT)*(xSC-xFT));
2506 /* FE 2019-10-31 in case separation could not be measured, i.e. 0, then use default scale */
2507 double pscale = scale;
2508 if (sts_mode != DUAL_STS) pscale = sep ? rho_in/sep : scale;
2509 /* double pscale = sep ? rho_in/sep : 0.; */
2510 qc_val=pscale;
2511 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2512 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2513 cpl_propertylist_set_comment (o_header, qc_name,
2514 "[mas/pixel] plate-scale in the "
2515 "FT-SC direction");
2516 }
2517
2518 /*------------------------------------------------*/
2519 /* If in dual field, measure fibre position error */
2520 /*------------------------------------------------*/
2521
2522 if (axis_mode != MODE_ONAXIS) {
2523 /* Error in SC fibre positioning */
2524 /* The three terms are */
2525 /* - offset from FT target as detected to original SC */
2526 /* target as detected; */
2527 /* - blind offset command from mapping template projected */
2528 /* on acqisition camera; */
2529 /* - offset from FT fiber to SC fiber. */
2530 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DX", tel+1);
2531 qc_val = 0;
2532 if (scale) {
2533 qc_val = (xSC-xFT) + sobj_offx_cam/scale - fiber_ft_sc_x - calib_dx[tel];
2534 }
2535 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2536 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2537 cpl_propertylist_set_comment (o_header, qc_name,
2538 "[pixel] dx from SC fiber to "
2539 "SC object");
2540
2541 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DY", tel+1);
2542 qc_val = 0;
2543 if (scale) {
2544 qc_val = (ySC-yFT) + sobj_offy_cam/scale - fiber_ft_sc_y - calib_dy[tel];
2545 }
2546 cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
2547 cpl_propertylist_update_double (o_header, qc_name, qc_val);
2548 cpl_propertylist_set_comment (o_header, qc_name,
2549 "[pixel] dx from SC fiber to "
2550 "SC object");
2551 }
2552
2553 /*-----------------------------*/
2554 /* Average Strehl computation */
2555 /*-----------------------------*/
2556
2557 double max_on_average, strehl_on_average;
2558
2559 /* Measure strehl and maximum on averaged image */
2560 gravi_acq_measure_strehl(mean_img, xFT, yFT, scale, &strehl_on_average, header);
2561
2562 /*
2563 * JIRA PIPE-9123 : Pipeline faile with Strehl=NaN
2564 * For 2017 AT data the spots jump up by 30 pixel and therefore are not
2565 * catched in the subwindow around roof_pos , cut_out +/- 50 pixel
2566 */
2567 if (cpl_error_get_code() != CPL_ERROR_NONE ) {
2568 cpl_msg_info (cpl_func, "WARNING Filling STREHL due to NaN");
2569 cpl_error_reset();
2570 strehl_on_average = 0.0;
2571 }
2572
2573 gravi_acq_measure_max(mean_img, xFT, yFT, 15, &max_on_average);
2574
2575 /* Update Strehl QC */
2576 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL", tel+1);
2577 cpl_msg_info (cpl_func, "%s = %f", qc_name, strehl_on_average);
2578 cpl_propertylist_update_double (o_header, qc_name, strehl_on_average);
2579 cpl_propertylist_set_comment (o_header, qc_name, "Average Strehl from stacked AcqCam images");
2580
2581
2582 /* Adding Strehl for SC channel in wide mode (PIPE-9913) */
2583 if (sts_mode == DUAL_STS) {
2584 gravi_acq_measure_strehl(mean_img, xSC, ySC, scale, &strehl_on_average, header);
2585
2586 /*
2587 * JIRA PIPE-9123 : Pipeline faile with Strehl=NaN
2588 * For 2017 AT data the spots jump up by 30 pixel and therefore are not
2589 * catched in the subwindow around roof_pos , cut_out +/- 50 pixel
2590 */
2591 if (cpl_error_get_code() != CPL_ERROR_NONE ) {
2592 cpl_msg_info (cpl_func, "WARNING Filling STREHL SC due to NaN");
2593 cpl_error_reset();
2594 strehl_on_average = 0.0;
2595 }
2596
2597 /* Update Strehl QC */
2598 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHLSC", tel+1);
2599 cpl_msg_info (cpl_func, "%s = %f", qc_name, strehl_on_average);
2600 cpl_propertylist_update_double (o_header, qc_name, strehl_on_average);
2601 cpl_propertylist_set_comment (o_header, qc_name, "Average Strehl from stacked AcqCam images");
2602 }
2603 /*----------------------------*/
2604 /* Now process frame by frame */
2605 /*----------------------------*/
2606
2607 /* Loop on all images */
2608 for (cpl_size row = 0; row < nrow; row++) {
2609 if (row %10 == 0 || row == (nrow-1))
2610 cpl_msg_info_overwritable (cpl_func, "Fit image %lld over %lld", row+1, nrow);
2611
2612 /* Get data */
2613 cpl_image * img = cpl_imagelist_get (acqcam_imglist, row);
2614 CPLCHECK_MSG("Error getting image");
2615
2616 /*-------------------------------------------------*/
2617 /* SC target position computation of current frame */
2618 /*-------------------------------------------------*/
2619 double xsc = xSC, ysc = ySC, exsc=0., eysc=0.;
2620 gravi_acq_fit_gaussian (img, &xsc, &ysc, &exsc, &eysc, size);
2621 CPLCHECK_MSG("Error fitting SC");
2622 /* Shift back positions to full frame */
2623 if (xsc != 0.) xsc += sx - 1 - nsx*tel;
2624 if (ysc != 0.) ysc += sy - 1;
2625 cpl_table_set (acqcam_table, "FIELD_SC_X", row*ntel+tel, xsc);
2626 cpl_table_set (acqcam_table, "FIELD_SC_Y", row*ntel+tel, ysc);
2627 cpl_table_set (acqcam_table, "FIELD_SC_XERR", row*ntel+tel, exsc);
2628 cpl_table_set (acqcam_table, "FIELD_SC_YERR", row*ntel+tel, eysc);
2629 CPLCHECK_MSG("Error setting SC columns");
2630
2631 /*-------------------------------------------------*/
2632 /* FT target position computation of current frame */
2633 /*-------------------------------------------------*/
2634 double xft = xFT, yft = yFT, exft=0., eyft=0.;
2635 if (axis_mode != MODE_ONAXIS) {
2636 gravi_acq_fit_gaussian (img, &xft, &yft, &exft, &eyft, size);
2637 CPLCHECK_MSG("Error fitting FT");
2638 /* Shift back positions to full frame */
2639 if (xft != 0.) xft += sx - 1 - nsx*tel;
2640 if (yft != 0.) yft += sy - 1;
2641 } else {
2642 xft=xsc;
2643 yft=ysc;
2644 exft=exsc;
2645 eyft=eysc;
2646 }
2647 cpl_table_set (acqcam_table, "FIELD_FT_X", row*ntel+tel, xft);
2648 cpl_table_set (acqcam_table, "FIELD_FT_Y", row*ntel+tel, yft);
2649 cpl_table_set (acqcam_table, "FIELD_FT_XERR", row*ntel+tel, exft);
2650 cpl_table_set (acqcam_table, "FIELD_FT_YERR", row*ntel+tel, eyft);
2651 CPLCHECK_MSG("Error setting FT column");
2652
2653 /*------------------------------------------*/
2654 /* Plate scale computation of current frame */
2655 /*------------------------------------------*/
2656 double ft_sc_x = xsc - xft;
2657 double ft_sc_y = ysc - yft;
2658 double eft_sc_x = sqrt(exsc*exsc+exft*exft);
2659 double eft_sc_y = sqrt(eysc*eysc+eyft*eyft);
2660 double sep = sqrt(ft_sc_x*ft_sc_x+ft_sc_y*ft_sc_y);
2661 double pscale = scale;
2662 if (sts_mode != DUAL_STS) pscale = sep ? rho_in/sep : 0.;
2663 double escale = 0.;
2664 if (sep) {
2665 escale = rho_in/(sep*sep*sep)*
2666 (ft_sc_x*eft_sc_x+ft_sc_y*eft_sc_y) ;
2667 }
2668 cpl_table_set (acqcam_table, "FIELD_SCALE", row*ntel+tel, pscale);
2669 cpl_table_set (acqcam_table, "FIELD_SCALEERR", row*ntel+tel, escale);
2670
2671 /*---------------------------------------------*/
2672 /* Error in fibre positioning of current frame */
2673 /*---------------------------------------------*/
2674 /* The three terms are */
2675 /* - offset from FT target as detected to original SC */
2676 /* target as detected; */
2677 /* - blind offset command from mapping template projected */
2678 /* on acqisition camera; */
2679 /* - offset from FT fiber to SC fiber. */
2680 double corrx=0, corry=0., ecorrx=0., ecorry=0.;
2681 if (pscale) {
2682 corrx = ft_sc_x + sobj_offx_cam/pscale - fiber_ft_sc_x - calib_dx[tel];
2683 corry = ft_sc_y + sobj_offy_cam/pscale - fiber_ft_sc_y - calib_dy[tel];
2684 double tmp = escale/(pscale*pscale);
2685 tmp *= tmp;
2686 ecorrx = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offx_cam*sobj_offx_cam*tmp);
2687 ecorry = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offy_cam*sobj_offy_cam*tmp);
2688 }
2689 cpl_table_set (acqcam_table, "FIELD_FIBER_DX", row*ntel+tel, corrx);
2690 cpl_table_set (acqcam_table, "FIELD_FIBER_DY", row*ntel+tel, corry);
2691 cpl_table_set (acqcam_table, "FIELD_FIBER_DXERR", row*ntel+tel, ecorrx);
2692 cpl_table_set (acqcam_table, "FIELD_FIBER_DYERR", row*ntel+tel, ecorry);
2693
2694 /*-------------------------------------*/
2695 /* Strehl computation of current frame */
2696 /*-------------------------------------*/
2697 double max_on_frame;
2698 gravi_acq_measure_max(img, xFT, yFT, 15, &max_on_frame);
2699 cpl_table_set (acqcam_table, "FIELD_STREHL", row*ntel+tel, strehl_on_average*(max_on_frame/max_on_average) );
2700
2701 } /* End loop on images */
2702
2703 /* Add some QC */
2704 double sc_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_SC_X", tel, ntel);
2705 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X STD", tel+1);
2706 cpl_propertylist_update_double (o_header, qc_name, sc_std_x);
2707 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
2708
2709 double sc_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_SC_Y", tel, ntel);
2710 sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y STD", tel+1);
2711 cpl_propertylist_update_double (o_header, qc_name, sc_std_y);
2712 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
2713
2714 double ft_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_FT_X", tel, ntel);
2715 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X STD", tel+1);
2716 cpl_propertylist_update_double (o_header, qc_name, ft_std_x);
2717 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
2718
2719 double ft_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_FT_Y", tel, ntel);
2720 sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y STD", tel+1);
2721 cpl_propertylist_update_double (o_header, qc_name, ft_std_y);
2722 cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
2723
2724 double strehl_std = gravi_table_get_column_std (acqcam_table, "FIELD_STREHL", tel, ntel);
2725 sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL STD", tel+1);
2726 cpl_propertylist_update_double (o_header, qc_name, strehl_std);
2727 cpl_propertylist_set_comment (o_header, qc_name, "Std of FT strehl");
2728
2729 } /* End loop on tel */
2730
2732 return CPL_ERROR_NONE;
2733}
2734
2735/*----------------------------------------------------------------------------*/
2752/*----------------------------------------------------------------------------*/
2753
2754cpl_error_code gravi_reduce_acqcam (gravi_data * output_data,
2755 gravi_data * input_data,
2756 gravi_data * sky_data,
2757 gravi_data * dark_data,
2758 gravi_data * static_param_data)
2759{
2761 cpl_ensure_code (output_data, CPL_ERROR_NULL_INPUT);
2762 cpl_ensure_code (input_data, CPL_ERROR_NULL_INPUT);
2763
2764 /* Check if extension exist */
2766 gravi_msg_warning (cpl_func, "Cannot reduce the ACQCAM, not data");
2767 return CPL_ERROR_NONE;
2768 }
2769
2770 /* Get the data and header */
2771 cpl_propertylist * header, * o_header;
2772 header = gravi_data_get_header (input_data);
2773 o_header = gravi_data_get_header (output_data);
2774
2775 /* Check if ISS data present (to avoid crash) */
2776 if (!gravi_conf_get_telname (0, header)) {
2777 gravi_msg_warning (cpl_func, "Cannot reduce the ACQCAM, no ISS keywords");
2778 return CPL_ERROR_NONE;
2779 }
2780
2781 cpl_imagelist * acqcam_imglist;
2782 cpl_imagelist * acqcam_imglist_v2;
2783 acqcam_imglist = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2784 acqcam_imglist_v2 = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2785 CPLCHECK_MSG ("Cannot get data or header");
2786
2787 /* Build the table */
2788 cpl_size ntel = 4;
2789 cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
2790
2791 cpl_table * acqcam_table;
2792 acqcam_table = cpl_table_new (nrow * ntel);
2793
2794 /* Time column shall contain the time from PCR.ACQ.START in [us] */
2795 cpl_table_new_column (acqcam_table, "TIME", CPL_TYPE_INT);
2796 cpl_table_set_column_unit (acqcam_table, "TIME", "us");
2797
2798 /* Loop on DIT in cube to fill the TIME column
2799 * same value for all 4 beams*/
2800 for (cpl_size row = 0; row < nrow; row++) {
2801 double time = gravi_pfits_get_time_acqcam (header, row);
2802 for (int tel = 0; tel < ntel; tel ++)
2803 cpl_table_set (acqcam_table, "TIME", row*ntel+tel, time);
2804 }
2805
2806 /* Compute mean image */
2807 cpl_image * mean_img = cpl_imagelist_collapse_create (acqcam_imglist);
2808 cpl_image * mean_subtracted_img = NULL;
2809
2810 if (sky_data) {
2811 mean_subtracted_img = gravi_data_get_img(sky_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2812 cpl_image_subtract(mean_subtracted_img, mean_img);
2813 cpl_image_multiply_scalar(mean_subtracted_img, -1.0);
2814 } else if (dark_data) {
2815 mean_subtracted_img = gravi_data_get_img(dark_data, GRAVI_IMAGING_DATA_ACQ_EXT);
2816 cpl_msg_warning(cpl_func, "No SKY, using DARK for acquisition field subtraction");
2817 cpl_image_subtract(mean_subtracted_img, mean_img);
2818 cpl_image_multiply_scalar(mean_subtracted_img, -1.0);
2819 } else {
2820 cpl_msg_warning(cpl_func, "No SKY or DARK, skipping acquisition field subtraction");
2821 mean_subtracted_img = mean_img;
2822 }
2823
2824 /* Compute FIELD columns */
2825 gravi_acqcam_field (mean_subtracted_img, acqcam_imglist, header,
2826 acqcam_table, o_header, static_param_data);
2827 CPLCHECK_MSG ("Cannot reduce acquisition field images");
2828
2829 /* Compute PUPIL columns with algorithm 2.0*/
2830 gravi_acqcam_pupil_v2 (mean_img, acqcam_imglist_v2, header,
2831 acqcam_table, o_header, static_param_data);
2832 CPLCHECK_MSG ("Cannot reduce pupil images");
2833
2834 /* Add this output table in the gravi_data */
2835 cpl_propertylist * plist_acq_cam = cpl_propertylist_new ();
2836 cpl_propertylist_update_string (plist_acq_cam, "INSNAME", INSNAME_ACQ);
2837 gravi_data_add_img (output_data, plist_acq_cam, GRAVI_IMAGING_DATA_ACQ_EXT, mean_img);
2838 CPLCHECK_MSG ("Cannot add acqcam_table in data");
2839
2840 gravi_data_add_table (output_data, NULL, GRAVI_OI_VIS_ACQ_EXT, acqcam_table);
2841 CPLCHECK_MSG ("Cannot add acqcam_table in data");
2842
2844 return CPL_ERROR_NONE;
2845}
2846
2847/*----------------------------------------------------------------------------*/
2852/*----------------------------------------------------------------------------*/
2853
2854double gravi_acqcam_z2meter (double PositionPixels, gravi_data * static_param_data)
2855{
2856 double f_PT ; /* = 14e-3; */ /* pupil tracker lenslet FL*/
2857 double f_lens ; /* = 467e-3; */ /* folding optics lens FL */
2858 double Llambda ; /* = 1.2e-6; */ /* laser diode wavelength */
2859 double D_beam ; /* = 18e-3; */ /* meter */
2860 double D_pixel ; /* = 18e-6; */
2861 double D_AT ; /* = 1.8; */ /* m */
2862 double D_lenslet ; /* = 2 * 1.015e-3; */
2863
2864 f_PT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS fPT");
2865 cpl_msg_debug (cpl_func,"ESO LENS fPT is : %e", f_PT);
2866
2867 f_lens = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS flens");
2868 cpl_msg_debug (cpl_func,"ESO LENS f_lens is : %e", f_lens);
2869
2870 Llambda = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Llambda");
2871 cpl_msg_debug (cpl_func,"ESO LENS Llambda is : %e", Llambda);
2872
2873 D_beam = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dbeam");
2874 cpl_msg_debug (cpl_func,"ESO LENS D_beam is : %e", D_beam);
2875
2876 D_pixel = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dpixel");
2877 cpl_msg_debug (cpl_func,"ESO LENS D_pixel is : %e", D_pixel);
2878
2879 D_AT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS DAT");
2880 cpl_msg_debug (cpl_func,"ESO LENS D_AT is : %e", D_AT);
2881
2882 D_lenslet = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dlenslet");
2883 cpl_msg_debug (cpl_func,"ESO LENS D_lenslet is : %e", D_lenslet);
2884
2885
2886 double longDef;
2887 longDef = 8 * (f_PT / D_lenslet) * (f_PT / D_lenslet) * 3.5 * D_pixel *
2888 D_beam / (f_PT * D_lenslet) * Llambda / CPL_MATH_2PI * PositionPixels;
2889
2890 return f_lens * f_lens * longDef / (f_PT + longDef) / f_PT * (D_AT / D_lenslet);
2891}
2892
2893
2894/*----------------------------------------------------------------------------*/
2905/*----------------------------------------------------------------------------*/
2906
2907cpl_error_code gravi_image_fft_correlate (cpl_image *ia, cpl_image *ib, cpl_size *xd, cpl_size *yd)
2908{
2910 cpl_ensure_code (ia, CPL_ERROR_NULL_INPUT);
2911 cpl_ensure_code (ib, CPL_ERROR_NULL_INPUT);
2912 cpl_ensure_code (xd, CPL_ERROR_NULL_INPUT);
2913 cpl_ensure_code (yd, CPL_ERROR_NULL_INPUT);
2914
2915 cpl_error_code error = CPL_ERROR_NONE;
2916 cpl_type type = cpl_image_get_type(ia);
2917 cpl_size nx = cpl_image_get_size_x(ia);
2918 cpl_size ny = cpl_image_get_size_y(ia);
2919
2920 cpl_image * ic = cpl_image_new(nx, ny, type);
2921 cpl_image * fa = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2922 cpl_image * fb = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2923 cpl_image * fc = cpl_image_new(nx, ny, type | CPL_TYPE_COMPLEX);
2924
2925 cpl_imagelist * iab = cpl_imagelist_new();
2926 cpl_imagelist * fab = cpl_imagelist_new();
2927
2928 /* Put image ia into imagelist iab */
2929 error = cpl_imagelist_set(iab, ia, 0);
2930
2931 /* Put image ib into imagelist iab */
2932 error = cpl_imagelist_set(iab, ib, 1);
2933
2934 /* Put empty Fourier image fa in imagelist fab */
2935 error = cpl_imagelist_set(fab, fa, 0);
2936
2937 /* Put empty Fourier image fb in imagelist fab*/
2938 error = cpl_imagelist_set(fab, fb, 1);
2939
2940 /* FFT on the image list iab and fill FFT image list fab */
2941 error = cpl_fft_imagelist(fab, iab, CPL_FFT_FORWARD);
2942
2943 /* conjugate fourier image list fb to fc */
2944 error = cpl_image_conjugate(fc, fb);
2945
2946 /* Multiply fourier image list fa and fc*/
2947 error = cpl_image_multiply(fc, fa);
2948
2949 /* Reverse FFT */
2950 error = cpl_fft_image(ic, fc, CPL_FFT_BACKWARD | CPL_FFT_NOSCALE);
2951
2952 /* Find max position -> shift*/
2953 error = cpl_image_get_maxpos(ic, xd, yd);
2954
2955 /* Unwrap for "negative" maximum position */
2956 if (*xd > nx/2) *xd = *xd - nx;
2957 if (*yd > ny/2) *yd = *yd - ny;
2958
2959 /* subtract one for come from fits convention (lower left pixel = 1,1) to shift vector */
2960 (*xd)--;
2961 (*yd)--;
2962
2963 // cpl_image_save (ia, "ia.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2964 // cpl_image_save (ib, "ib.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2965 // cpl_image_save (ic, "ic.fits", CPL_TYPE_DOUBLE,NULL, CPL_IO_DEFAULT);
2966
2967 /* Free memory allocated by this routine */
2968 FREE (cpl_imagelist_unwrap, iab);
2969 FREE (cpl_imagelist_unwrap, fab);
2970 FREE (cpl_image_delete, ic);
2971 FREE (cpl_image_delete, fa);
2972 FREE (cpl_image_delete, fb);
2973 FREE (cpl_image_delete, fc);
2974
2976 return (error);
2977}
2978
#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_reduce_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *sky_data, gravi_data *dark_data, gravi_data *static_param_data)
Reduce the ACQ camera images.
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:946
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:817
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_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
Definition: gravi_cpl.c:1656
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:2004
double gravi_table_get_column_std(cpl_table *table, const char *name, int base, int nbase)
Definition: gravi_cpl.c:388
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