GRAVI Pipeline Reference Manual  1.2.3
gravi_acqcam.c
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  */
41 /*-----------------------------------------------------------------------------
42  Includes
43  -----------------------------------------------------------------------------*/
44 
45 #ifdef HAVE_CONFIG_H
46 #include <config.h>
47 #endif
48 
49 #include <cpl.h>
50 #include <string.h>
51 #include <stdio.h>
52 #include <math.h>
53 #include <time.h>
54 
55 #include <hdrl_strehl.h>
56 
57 #include "gravi_dfs.h"
58 #include "gravi_data.h"
59 #include "gravi_pfits.h"
60 
61 #include "gravi_cpl.h"
62 #include "gravi_utils.h"
63 
64 #include "gravi_acqcam.h"
65 
66 /*-----------------------------------------------------------------------------
67  Private prototypes
68  -----------------------------------------------------------------------------*/
69 
70 double exp1 (double x);
71 double sin1 (double x);
72 int gravi_acqcam_xy_diode (const double v[], double *xd, double *yd);
73 
74 cpl_error_code gravi_acqcam_get_pup_ref (cpl_propertylist * header,
75  cpl_size tel,
76  cpl_vector * pupref);
77 
78 cpl_error_code gravi_acqcam_get_diode_ref (cpl_propertylist * header,
79  cpl_size tel,
80  cpl_vector * output);
81 
82 static int gravi_acqcam_spot (const double x_in[], const double v[], double *result);
83 static int gravi_acqcam_xy_sub (const double v[], double *xsub, double *ysub);
84 
85 cpl_error_code gravi_acqcam_spot_imprint (cpl_image * img, cpl_vector * a);
86 
87 static int gravi_acqcam_spot_dfda (const double x_in[], const double v[], double result[]);
88 
89 cpl_error_code gravi_acqcam_fit_spot (cpl_image * img, cpl_size ntry,
90  cpl_vector * a,
91  int fitAll,
92  int * nspot);
93 
94 double gravi_acqcam_z2meter (double PositionPixels, gravi_data *static_param_data);
95 
96 cpl_error_code gravi_acqcam_pupil (cpl_image * mean_img,
97  cpl_imagelist * acqcam_imglist,
98  cpl_propertylist * header,
99  cpl_table * acqcam_table,
100  cpl_propertylist * o_header,
101  gravi_data *static_param_data);
102 
103 cpl_error_code gravi_acqcam_field (cpl_image * mean_img,
104  cpl_imagelist * acqcam_imglist,
105  cpl_propertylist * header,
106  cpl_table * acqcam_table,
107  cpl_propertylist * o_header,
108  gravi_data *static_param_data);
109 
110 cpl_error_code gravi_acq_fit_gaussian (cpl_image * img, double *x, double *y,
111  double *ex, double *ey, cpl_size size);
112 
113 cpl_error_code gravi_acq_measure_strehl(cpl_image * img, double x, double y,
114  double pscale, double *SR, cpl_propertylist * header);
115 
116 cpl_error_code gravi_acq_measure_max(cpl_image * img, double x, double y, double size, double * img_max);
117 
118 /* This global variable optimises the computation
119  * of partial derivative on fitted parameters */
120 const extern int * GRAVI_LVMQ_FREE;
121 const int * GRAVI_LVMQ_FREE = NULL;
122 
123 /* Number of parameters in the model 'gravi_acqcam_spot'
124  * And position of parameters */
125 #define GRAVI_SPOT_NA 30
126 #define GRAVI_SPOT_SUB 0
127 #define GRAVI_SPOT_ANGLE 8
128 #define GRAVI_SPOT_SCALE 9
129 #define GRAVI_SPOT_DIODE 10
130 #define GRAVI_SPOT_FWHM 13
131 #define GRAVI_SPOT_FLUX 14
132 
133 #define GRAVI_ACQ_PUP_FLUX 1e6
134 
135 /*-----------------------------------------------------------------------------
136  Functions code
137  -----------------------------------------------------------------------------*/
138 
139 /*----------------------------------------------------------------------------*/
152 /*----------------------------------------------------------------------------*/
153 
154 cpl_error_code gravi_preproc_acqcam (gravi_data *output_data,
155  gravi_data *input_data,
156  gravi_data *bad_map)
157 {
158  gravi_msg_function_start(1);
159  cpl_ensure_code (output_data, CPL_ERROR_NULL_INPUT);
160  cpl_ensure_code (input_data, CPL_ERROR_NULL_INPUT);
161 
162  /* Check if extension exist */
163  if (!gravi_data_has_extension (input_data, GRAVI_IMAGING_DATA_ACQ_EXT)) {
164  gravi_msg_warning (cpl_func,"Cannot preproc the ACQCAM, not data in file");
165  return CPL_ERROR_NONE;
166  }
167  if (!gravi_data_has_extension (bad_map, GRAVI_IMAGING_DATA_ACQ_EXT)) {
168  gravi_msg_warning (cpl_func,"Cannot preproc the ACQCAM, no badpixel in BAD");
169  return CPL_ERROR_NONE;
170  }
171 
172  /* Construct a mask of badpixels */
173  cpl_image * badpix_img = gravi_data_get_img (bad_map, GRAVI_IMAGING_DATA_ACQ_EXT);
174  cpl_mask * badpix_mask = cpl_mask_threshold_image_create (badpix_img, 0.5, 10000);
175  CPLCHECK_MSG ("Cannot get BAD map for ACQ");
176 
177  /* Get the imagelist */
178  cpl_imagelist * imglist;
179  imglist = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
180  CPLCHECK_MSG ("Cannot get image for ACQ");
181 
182  /* Allocate new memory */
183  imglist = cpl_imagelist_duplicate (imglist);
184 
185  /*
186  * Loop on images to cleanup-badpixels
187  */
188  cpl_size nrow = cpl_imagelist_get_size (imglist);
189  for (cpl_size row = 0; row < nrow; row++) {
190 
191  /* Get image */
192  cpl_image * img = cpl_imagelist_get (imglist, row);
193 
194  /* Cleanup-badpixel */
195  cpl_image_reject_from_mask (img, badpix_mask);
196  cpl_detector_interpolate_rejected (img);
197  CPLCHECK_MSG ("Cannot clean-up badpixel");
198  }
199 
200  FREE (cpl_mask_delete, badpix_mask);
201 
202  /* Get the size */
203  cpl_image * img = cpl_imagelist_get (imglist, 0);
204  cpl_size nx = cpl_image_get_size_x (img);
205  cpl_size ny = cpl_image_get_size_y (img);
206 
207  /* FIXME: Deal with the case full-frame in a better way */
208  cpl_size ury = (ny>1100) ? 1200 : 750;
209 
210  /*
211  * Remove the pupil background by the mean of blinking
212  */
213 
214  if (nrow == 1) {
215  gravi_msg_warning ("FIXME","Cannot remove blinked pupil (no blink)");
216 
217  } else {
218  cpl_msg_info (cpl_func, "Remove the blinking");
219 
220  /* Get pupil diode flux for all images. Deal with the
221  * change of format (FIXME: make it more clean) */
222  double flux[nrow];
223  for (cpl_size row = 0; row < nrow; row ++) {
224  cpl_image * img = cpl_imagelist_get (imglist,row);
225  if (nx < 1100)
226  flux[row] = cpl_image_get_flux_window (img, 1, 800, 1000, 1000) -
227  cpl_image_get_median_window (img, 1, 800, 1000, 1000) * 1000 * 200;
228  else
229  flux[row] = cpl_image_get_flux_window (img, 1, 1200, 2048, 1536) -
230  cpl_image_get_median_window (img, 1, 1200, 2048, 1536) * 2048 * 336;
231 
232  cpl_msg_debug (cpl_func, "flux %lli = %.2f",row, flux[row]);
233  }
234 
235  for (cpl_size row = 0; row < nrow; row ++) {
236 
237  /* Find the best possible blink, default is
238  * current (hence frame will be zero) */
239  cpl_size blink = row;
240  if ( (flux[row] - flux[CPL_MAX(row-1,0)] ) > GRAVI_ACQ_PUP_FLUX) blink = row-1;
241  if ( (flux[row] - flux[CPL_MIN(row+1,nrow-1)] ) > GRAVI_ACQ_PUP_FLUX) blink = row+1;
242 
243  cpl_msg_debug (cpl_func, "row %lli is debiased with %lli",row,blink);
244 
245  /* Remove the blink only for pupil part of the camera */
246  gravi_image_subtract_window (cpl_imagelist_get (imglist, row),
247  cpl_imagelist_get (imglist, blink),
248  1, ury, nx, ny, 1, ury);
249  CPLCHECK_MSG ("Cannot remove blinked pupil");
250  }
251  }
252 
253  /*
254  * Run a median filter for bias
255  */
256  cpl_msg_info (cpl_func, "Remove a median filter");
257 
258  cpl_mask * kernel = cpl_mask_new (11, 11);
259  cpl_mask_not (kernel);
260 
261  for (cpl_size row = 0; row < nrow; row ++) {
262  if (row %10 == 0)
263  cpl_msg_info_overwritable (cpl_func, "Median filter of ACQ %lld over %lld", row+1, nrow);
264 
265  cpl_image * img = cpl_imagelist_get (imglist, row);
266  cpl_image * unfiltered_img = cpl_image_extract (img, 1, ury, nx, ny);
267  cpl_image * filtered_img = cpl_image_duplicate (unfiltered_img);
268 
269  cpl_image_filter_mask (filtered_img, unfiltered_img, kernel,
270  CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
271  gravi_image_subtract_window (img, filtered_img,
272  1, ury, nx, ny, 1, 1);
273 
274  FREE (cpl_image_delete, unfiltered_img);
275  FREE (cpl_image_delete, filtered_img);
276  CPLCHECK_MSG ("Cannot run median filter");
277  }
278 
279  FREE (cpl_mask_delete, kernel);
280 
281  /*
282  * Set in output
283  */
284  gravi_data_add_cube (output_data, NULL, GRAVI_IMAGING_DATA_ACQ_EXT, imglist);
285 
286  gravi_msg_function_exit(1);
287  return CPL_ERROR_NONE;
288 }
289 
290 /*----------------------------------------------------------------------------*/
291 
292 /* Fast sin function */
293 double sin1 (double x)
294 {
295  /* Within -pi +pi*/
296  while (x >= CPL_MATH_PI) x -= CPL_MATH_2PI;
297  while (x <= -CPL_MATH_PI) x += CPL_MATH_2PI;
298 
299  /* Within -pi/2 +pi/2*/
300  if (x > CPL_MATH_PI_2) x = CPL_MATH_PI - x;
301  if (x < -CPL_MATH_PI_2) x = -CPL_MATH_PI - x;
302 
303  double x2 = x*x, x3 = x2*x, x5 = x3*x2;
304  return 0.9996949 * x - 0.1656700 * x3 + 0.0075134 * x5;
305 }
306 
307 /* Fast exp function */
308 double exp1 (double x) {
309  x = 1.0 + x / 256.0;
310  x *= x; x *= x; x *= x; x *= x;
311  x *= x; x *= x; x *= x; x *= x;
312  return x;
313 }
314 
315 /* Compute the 4 diode positions from {angle, scaling, dx, dy}
316  * The model assume the 4 diodes form a rectangle centered on
317  * the pupil. Hence this model has a 180deg degeneracy */
318 int gravi_acqcam_xy_diode (const double v[], double *xd, double *yd)
319 {
320  double dx = v[GRAVI_SPOT_DIODE+0];
321  double dy = v[GRAVI_SPOT_DIODE+1];
322 
323  /* Angle */
324  double ang = (v[GRAVI_SPOT_ANGLE] - v[GRAVI_SPOT_DIODE+2])* CPL_MATH_RAD_DEG;
325  double sang = sin1 (ang) * v[GRAVI_SPOT_SCALE];
326  double cang = sin1 (ang + CPL_MATH_PI_2) * v[GRAVI_SPOT_SCALE];
327 
328  /* Diode arrangement */
329  xd[0] = -cang * dx + sang * dy;
330  xd[1] = -cang * dx - sang * dy;
331  xd[2] = cang * dx - sang * dy;
332  xd[3] = cang * dx + sang * dy;
333 
334  yd[0] = cang * dy + sang * dx;
335  yd[1] = -cang * dy + sang * dx;
336  yd[2] = -cang * dy - sang * dx;
337  yd[3] = cang * dy - sang * dx;
338 
339  return 0;
340 }
341 
342 /* Compute the 4 sub-aperture positions from the sub-aperture modes */
343 static int gravi_acqcam_xy_sub (const double v[], double *xsub, double *ysub)
344 {
345  /* Sub-apperture arrangement */
346  const double * vd = v + GRAVI_SPOT_SUB;
347  xsub[0] = vd[0] + vd[1] + vd[2] + vd[3];
348  xsub[1] = vd[0] - vd[1] + vd[2] - vd[3];
349  xsub[2] = vd[0] + vd[1] - vd[2] - vd[3];
350  xsub[3] = vd[0] - vd[1] - vd[2] + vd[3];
351 
352  ysub[0] = vd[4] + vd[5] + vd[6] + vd[7];
353  ysub[1] = vd[4] - vd[5] + vd[6] - vd[7];
354  ysub[2] = vd[4] + vd[5] - vd[6] - vd[7];
355  ysub[3] = vd[4] - vd[5] - vd[6] + vd[7];
356 
357  return 0;
358 }
359 
360 /*----------------------------------------------------------------------------*/
361 
362 static int gravi_acqcam_spot (const double x_in[], const double v[], double *result)
363 {
364  *result = 0.0;
365 
366  /* Static parameters */
367  double xd[4], yd[4], xsub[4], ysub[4];
368  gravi_acqcam_xy_diode (v, xd, yd);
369  gravi_acqcam_xy_sub (v, xsub, ysub);
370 
371  double fwhm2 = v[GRAVI_SPOT_FWHM];
372 
373  /* Loop on diode and appertures.
374  * The capture range is 2.FWHM */
375  for (int diode = 0; diode < 4 ; diode++) {
376  for (int sub = 0; sub < 4 ; sub++) {
377  double xf = (x_in[0] - xsub[sub] - xd[diode]);
378  double yf = (x_in[1] - ysub[sub] - yd[diode]);
379  double dist = (xf*xf + yf*yf) / fwhm2;
380  if (dist < 4.0) *result += v[GRAVI_SPOT_FLUX+sub*4+diode] * exp1 (-dist);
381  }
382  }
383 
384  return 0;
385 }
386 
387 /*----------------------------------------------------------------------------*/
388 
389 static int gravi_acqcam_spot_dfda (const double x_in[], const double v[], double result[])
390 {
391  double next = 0.0, here = 0.0, epsilon = 1e-8;
392 
393  double vlocal[GRAVI_SPOT_NA];
394  memcpy (vlocal, v, sizeof(double)*GRAVI_SPOT_NA);
395 
396  /* Compute value in-place */
397  gravi_acqcam_spot (x_in, vlocal, &here);
398 
399  /* Fill with zeros */
400  for (int a = 0; a < GRAVI_SPOT_NA; a++) result[a] = 0.0;
401 
402  /* Loop on parameters to compute finite differences
403  * FIXME: The analytical derivative may be faster, but
404  * wasn't true in first tests, thus keep these. */
405  for (int a = 0; a < GRAVI_SPOT_FLUX; a++) {
406  if (GRAVI_LVMQ_FREE[a] != 0) {
407  vlocal[a] += epsilon;
408  gravi_acqcam_spot (x_in, vlocal, &next);
409  vlocal[a] -= epsilon;
410  result[a] = (next - here) / epsilon;
411  }
412  }
413 
414  /* The intensities are trivial analytic derivative */
415  for (int a = GRAVI_SPOT_FLUX; a < GRAVI_SPOT_NA; a++) {
416  if (GRAVI_LVMQ_FREE[a] != 0) {
417  result[a] = here / v[a];
418  }
419  }
420 
421  return 0;
422 }
423 
424 /*----------------------------------------------------------------------------*/
431 /*----------------------------------------------------------------------------*/
432 
433 cpl_error_code gravi_acqcam_spot_imprint (cpl_image * img, cpl_vector * a)
434 {
435  gravi_msg_function_start(0);
436  cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
437  cpl_ensure_code (a, CPL_ERROR_NULL_INPUT);
438 
439  cpl_size nx = cpl_image_get_size_x (img);
440  cpl_size ny = cpl_image_get_size_y (img);
441 
442  /* Static parameters */
443  double xd[4], yd[4], xsub[4], ysub[4];
444  double * v = cpl_vector_get_data (a);
445  gravi_acqcam_xy_diode (v, xd, yd);
446  gravi_acqcam_xy_sub (v, xsub, ysub);
447 
448  /* Loop on diode and appertures */
449  for (int diode = 0; diode < 4 ; diode++) {
450  for (int sub = 0; sub < 4 ; sub++) {
451  cpl_size xf = roundl(xsub[sub] + xd[diode]) + 1;
452  cpl_size yf = roundl(ysub[sub] + yd[diode]) + 1;
453  if (xf < 2 || xf > nx-2 || yf < 2 || yf > ny-2) continue;
454  cpl_image_set (img, xf, yf, 0);
455  cpl_image_set (img, xf-1, yf, 0);
456  cpl_image_set (img, xf+1, yf, 0);
457  cpl_image_set (img, xf, yf+1, 0);
458  cpl_image_set (img, xf, yf-1, 0);
459  CPLCHECK ("Cannot imprint cross in image");
460  }
461  }
462 
463  gravi_msg_function_exit(0);
464  return CPL_ERROR_NONE;
465 }
466 
467 /*----------------------------------------------------------------------------*/
486 /*----------------------------------------------------------------------------*/
487 
488 cpl_error_code gravi_acqcam_get_diode_ref (cpl_propertylist * header,
489  cpl_size tel,
490  cpl_vector * output)
491 {
492  gravi_msg_function_start(0);
493  cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
494  cpl_ensure_code (tel>=0 && tel<4, CPL_ERROR_ILLEGAL_INPUT);
495  cpl_ensure_code (output, CPL_ERROR_ILLEGAL_INPUT);
496 
497  /* Get the telescope name and ID */
498  const char * telname = gravi_conf_get_telname (tel, header);
499 
500  /* Check telescope name */
501  if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
502  cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
503 
504  /* Hardcoded theoretical positions in mm */
505 
506  /* If UTs or ATs, select scaling, rotation, and spacing */
507  if (telname[0] == 'U') {
508  cpl_vector_set (output, GRAVI_SPOT_ANGLE, 0.0);
509  cpl_vector_set (output, GRAVI_SPOT_SCALE, 16.225);
510  cpl_vector_set (output, GRAVI_SPOT_DIODE+0, 0.363);
511  cpl_vector_set (output, GRAVI_SPOT_DIODE+1, 0.823);
512  cpl_vector_set (output, GRAVI_SPOT_DIODE+2, 0);
513  } else if (telname[0] == 'A') {
514  cpl_vector_set (output, GRAVI_SPOT_ANGLE, 0.0);
515  cpl_vector_set (output, GRAVI_SPOT_SCALE, 73.0154);
516  cpl_vector_set (output, GRAVI_SPOT_DIODE+0, 0.122);
517  cpl_vector_set (output, GRAVI_SPOT_DIODE+1, 0.158);
518  cpl_vector_set (output, GRAVI_SPOT_DIODE+2, 0);
519  } else
520  return cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
521  "Cannot get telescope name");
522 
523  gravi_msg_function_exit(0);
524  return CPL_ERROR_NONE;
525 }
526 
527 /*----------------------------------------------------------------------------*/
548 /*----------------------------------------------------------------------------*/
549 
550 cpl_error_code gravi_acqcam_get_pup_ref (cpl_propertylist * header,
551  cpl_size tel,
552  cpl_vector * output)
553 {
554  gravi_msg_function_start(0);
555  cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
556  cpl_ensure_code (tel>=0 && tel<4, CPL_ERROR_ILLEGAL_INPUT);
557  cpl_ensure_code (output, CPL_ERROR_ILLEGAL_INPUT);
558 
559  cpl_size ntel = 4, nsub = 4;
560  cpl_size nsx = 0, nsy = 0, sx = 0, sy = 0;
561 
562  /* If sub-windowing, we read the sub-window size and
563  * the sub-window start for pupil */
564  if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
565  char name[90];
566 
567  nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
568  sprintf (name, "ESO DET1 FRAM%lld STRX", 3*ntel + tel + 1);
569  sx = cpl_propertylist_get_int (header, name);
570 
571  nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY");
572  sprintf (name, "ESO DET1 FRAM%lld STRY", 3*ntel + tel + 1);
573  sy = cpl_propertylist_get_int (header, name);
574 
575  CPLCHECK_MSG ("Cannot get sub-windowing parameters");
576  }
577 
578  cpl_msg_debug (cpl_func,"sub-window pupil %lli sx= %lld sy = %lld", tel, sx, sy);
579 
580  /* Read the sub-apperture reference positions
581  * Converted to accound for sub-windowing
582  * In vector convention (start at 0,0) */
583  double xsub[4], ysub[4];
584  for (int sub = 0; sub < nsub ; sub++) {
585  double xv = gravi_pfits_get_ptfc_acqcam (header, sub*ntel + tel + 1);
586  double yv = gravi_pfits_get_ptfc_acqcam (header, sub*ntel + tel + 17);
587  xsub[sub] = xv - (sx - tel*nsx) - 1;
588  ysub[sub] = yv - (sy - 3*nsy) - 1;
589  cpl_msg_debug (cpl_func,"pupil %lli subC %i = %10.4f,%10.4f",
590  tel, sub, xsub[sub], ysub[sub]);
591  CPLCHECK_MSG ("Cannot get pupil reference position");
592  }
593 
594  /* All linear combination of sub-appertures center */
595  cpl_vector_set (output, GRAVI_SPOT_SUB+0, 0.25 * (xsub[0] + xsub[1] + xsub[2] + xsub[3]));
596  cpl_vector_set (output, GRAVI_SPOT_SUB+1, 0.25 * (xsub[0] - xsub[1] + xsub[2] - xsub[3]));
597  cpl_vector_set (output, GRAVI_SPOT_SUB+2, 0.25 * (xsub[0] + xsub[1] - xsub[2] - xsub[3]));
598  cpl_vector_set (output, GRAVI_SPOT_SUB+3, 0.25 * (xsub[0] - xsub[1] - xsub[2] + xsub[3]));
599 
600  cpl_vector_set (output, GRAVI_SPOT_SUB+4, 0.25 * (ysub[0] + ysub[1] + ysub[2] + ysub[3]));
601  cpl_vector_set (output, GRAVI_SPOT_SUB+5, 0.25 * (ysub[0] - ysub[1] + ysub[2] - ysub[3]));
602  cpl_vector_set (output, GRAVI_SPOT_SUB+6, 0.25 * (ysub[0] + ysub[1] - ysub[2] - ysub[3]));
603  cpl_vector_set (output, GRAVI_SPOT_SUB+7, 0.25 * (ysub[0] - ysub[1] - ysub[2] + ysub[3]));
604 
605  gravi_msg_function_exit(0);
606  return CPL_ERROR_NONE;
607 }
608 
609 /*----------------------------------------------------------------------------*/
629 /*----------------------------------------------------------------------------*/
630 
631 cpl_error_code gravi_acqcam_fit_spot (cpl_image * img,
632  cpl_size ntry,
633  cpl_vector * a,
634  int fitAll,
635  int * nspot)
636 {
637  gravi_msg_function_start(0);
638  cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
639  cpl_ensure_code (a, CPL_ERROR_NULL_INPUT);
640  cpl_ensure_code (ntry>0, CPL_ERROR_ILLEGAL_INPUT);
641  cpl_ensure_code (nspot, CPL_ERROR_NULL_INPUT);
642 
643  int nv = 0;
644  cpl_size x0 = cpl_vector_get (a, GRAVI_SPOT_SUB+0);
645  cpl_size y0 = cpl_vector_get (a, GRAVI_SPOT_SUB+4);
646  CPLCHECK_MSG ("Cannot get values valid");
647 
648  /* Compute RMS in the central region */
649  double RMS = gravi_image_get_noise_window (img, x0-25, y0-25, x0+25, y0+25);
650 
651  /* The image is surely empty */
652  if (RMS == 0) { *nspot = 0; return CPL_ERROR_NONE;}
653 
654  /*
655  * Coarse: fit with a re-bin image
656  */
657 
658  /* To lower the number of point, we extract a window of 175x175
659  * around the center of sub-apertures, rebin with 5x5 pixels */
660  cpl_size nw = 175, n_mean = 5;
661  cpl_size nc = nw/n_mean, nint = n_mean*n_mean;
662 
663  /* Allocate vector and matrix for the fit */
664  cpl_matrix * x_matrix = cpl_matrix_new (nc*nc, 2);
665  cpl_vector * y_vector = cpl_vector_new (nc*nc);
666  cpl_vector * sy_vector = cpl_vector_new (nc*nc);
667 
668  /* Loop on pixel in the considered window */
669  for (cpl_size x = 0; x < nc; x++) {
670  for (cpl_size y = 0; y < nc; y++) {
671  /* Average 3x3 pixels */
672  double z_mean = 0.0, x_mean = 0.0, y_mean = 0.0;
673  for (int i=-2; i<=2; i++) {
674  for (int j=-2; j<=2; j++) {
675  cpl_size x1 = x*n_mean+i+x0-(n_mean*nc)/2;
676  cpl_size y1 = y*n_mean+j+y0-(n_mean*nc)/2;
677  z_mean += cpl_image_get (img, x1+1, y1+1, &nv);
678  x_mean += x1; y_mean += y1;
679  }
680  }
681  /* Set in matrix and vectors */
682  cpl_matrix_set (x_matrix, x*nc+y, 0, x_mean/nint);
683  cpl_matrix_set (x_matrix, x*nc+y, 1, y_mean/nint);
684  cpl_vector_set (y_vector, x*nc+y, z_mean/nint);
685  cpl_vector_set (sy_vector, x*nc+y, RMS);
686  CPLCHECK_MSG ("Cannot fill matrix/vector");
687  }
688  } /* End loop on re-sampled pixels*/
689 
690  /* Output for global minimisation */
691  cpl_vector * a_start = cpl_vector_duplicate (a);
692  cpl_vector * a_tmp = cpl_vector_duplicate (a);
693  double chisq_final = 1e10;
694  srand(1);
695 
696  /* Set the fwhm to 6 to force a large capture range */
697  cpl_vector_set (a_start, GRAVI_SPOT_FWHM, 6.*6.);
698 
699  /* If needed, define a realistic value for the amplitude (for numerical stability) */
700  if (fitAll) for (int d=0;d<16;d++) cpl_vector_set (a_start, GRAVI_SPOT_FLUX+d, RMS);
701 
702  /* Fit sub-aperture mean position; and diode rotation */
703  const int ia_global[] = {1,0,0,0, 1,0,0,0, 1,0, 0,0,0, 0,
704  0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
705  GRAVI_LVMQ_FREE = ia_global;
706 
707  /* Loop on various starting points */
708  for (cpl_size try = 0; try < ntry; try++) {
709 
710  /* Move starting point in position (+-10pix) and angle (entire circle) */
711  cpl_vector_copy (a_tmp, a_start);
712  if (try > 0) {
713  cpl_vector_set (a_tmp, GRAVI_SPOT_SUB+0,
714  cpl_vector_get (a_tmp, GRAVI_SPOT_SUB+0) +
715  (rand()%20) - 10);
716  cpl_vector_set (a_tmp, GRAVI_SPOT_SUB+4,
717  cpl_vector_get (a_tmp, GRAVI_SPOT_SUB+4) +
718  (rand()%20) - 10);
719  cpl_vector_set (a_tmp, GRAVI_SPOT_ANGLE,
720  cpl_vector_get (a_tmp, GRAVI_SPOT_ANGLE) +
721  (rand()%180));
722  }
723 
724  /* Fit from this starting point */
725  double chisq;
726  cpl_fit_lvmq (x_matrix, NULL, y_vector, sy_vector,
727  a_tmp, GRAVI_LVMQ_FREE,
728  gravi_acqcam_spot,
729  gravi_acqcam_spot_dfda,
730  CPL_FIT_LVMQ_TOLERANCE,
731  CPL_FIT_LVMQ_COUNT,
732  CPL_FIT_LVMQ_MAXITER,
733  NULL, &chisq, NULL);
734  CPLCHECK_MSG ("Cannot fit global minimum");
735 
736  /* Check chi2 and keep if lowest so far */
737  if (chisq < chisq_final) {
738  cpl_msg_debug (cpl_func, "chisq_final = %f -> %f", chisq_final, chisq);
739  chisq_final = chisq;
740  cpl_vector_copy (a, a_tmp);
741  }
742 
743  } /* End loop on try starting points */
744 
745  FREE (cpl_matrix_delete, x_matrix);
746  FREE (cpl_vector_delete, y_vector);
747  FREE (cpl_vector_delete, sy_vector);
748  FREE (cpl_vector_delete, a_tmp);
749  FREE (cpl_vector_delete, a_start);
750 
751 
752  /*
753  * Fine: fit 10 pixel around each spot with true Gaussian
754  */
755 
756  cpl_size nf = 10, ndiode = 4, nsub = 4;
757  x_matrix = cpl_matrix_new (nf*nf*ndiode*nsub, 2);
758  y_vector = cpl_vector_new (nf*nf*ndiode*nsub);
759  sy_vector = cpl_vector_new (nf*nf*ndiode*nsub);
760 
761  double xd[4], yd[4], xsub[4], ysub[4];
762  gravi_acqcam_xy_diode (cpl_vector_get_data (a), xd, yd);
763  gravi_acqcam_xy_sub (cpl_vector_get_data (a), xsub, ysub);
764 
765  /* Loop on diode and appertures to fill the matrix and
766  * vector for the fine fit */
767  for (cpl_size v = 0, diode = 0; diode < ndiode ; diode++) {
768  for (int sub = 0; sub < nsub ; sub++) {
769  cpl_size xf = roundl(xsub[sub] + xd[diode]);
770  cpl_size yf = roundl(ysub[sub] + yd[diode]);
771 
772  /* Extract 10 pixels around each spot */
773  double mf = 0.0;
774  for (cpl_size x = xf-nf/2; x < xf+nf/2; x++) {
775  for (cpl_size y = yf-nf/2; y < yf+nf/2; y++) {
776  double value = cpl_image_get (img, x+1, y+1, &nv);
777  if (value > mf) mf = value;
778  cpl_matrix_set (x_matrix, v, 0, x);
779  cpl_matrix_set (x_matrix, v, 1, y);
780  cpl_vector_set (y_vector, v, value);
781  cpl_vector_set (sy_vector, v, RMS);
782  v++;
783  CPLCHECK_MSG ("Cannot fill matrix");
784  }
785  }
786 
787  /* Save the local maximum of this spot as
788  * the starting point for its amplitude */
789  cpl_vector_set (a, GRAVI_SPOT_FLUX+sub*4+diode, mf);
790  }
791  }
792 
793  /* Set FWHM to a realist value in [pix**2] */
794  cpl_vector_set (a, GRAVI_SPOT_FWHM, 2.3*2.3);
795 
796  /* Fit all sub-aperture position; rotation and scaling of diodes;
797  * and individual intensities of spots. (high order and scaling
798  * are only fitted if fitAll != 0) */
799  int F = fitAll;
800  const int ia_fine[] = {1,F,1,F, 1,1,F,F, 1,F, 0,0,0, 0,
801  1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
802  GRAVI_LVMQ_FREE = ia_fine;
803 
804  /* Fit from this starting point */
805  double chisq_fine = 0.0;
806  cpl_errorstate prestate = cpl_errorstate_get();
807  cpl_fit_lvmq (x_matrix, NULL, y_vector, sy_vector,
808  a, GRAVI_LVMQ_FREE,
809  gravi_acqcam_spot,
810  gravi_acqcam_spot_dfda,
811  CPL_FIT_LVMQ_TOLERANCE,
812  CPL_FIT_LVMQ_COUNT,
813  CPL_FIT_LVMQ_MAXITER,
814  NULL, &chisq_fine, NULL);
815 
816  /* Recover on error but return NULL */
817  if (cpl_error_get_code() == CPL_ERROR_CONTINUE) {
818  cpl_errorstate_set (prestate);
819  cpl_msg_warning (cpl_func, "Cannot fit pupil spots...");
820  for (int d=0; d<16; d++) cpl_vector_set (a, GRAVI_SPOT_FLUX+d, 0);
821  }
822 
823  CPLCHECK_MSG ("Cannot fit fine");
824 
825  cpl_msg_debug (cpl_func, "chisq_final = %.2f -> fine = %.2f",
826  chisq_final, chisq_fine);
827 
828  FREE (cpl_matrix_delete, x_matrix);
829  FREE (cpl_vector_delete, y_vector);
830  FREE (cpl_vector_delete, sy_vector);
831 
832  /* Count the number of valid spots based on their amplitude */
833  *nspot = 0;
834  for (int d=0; d<16; d++) if (cpl_vector_get (a, GRAVI_SPOT_FLUX+d) > 3.*RMS) (*nspot)++;
835 
836  gravi_msg_function_exit(0);
837  return CPL_ERROR_NONE;
838 }
839 
840 /*----------------------------------------------------------------------------*/
841 /*
842  * @brief measure Strehl Ratio of the source at the given location
843  */
844 /*----------------------------------------------------------------------------*/
845 
846 cpl_error_code gravi_acq_measure_strehl(cpl_image * img, double x, double y,
847  double pscale, double *SR, cpl_propertylist * header)
848 {
849  gravi_msg_function_start(0);
850  cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
851  cpl_ensure_code (SR, CPL_ERROR_NULL_INPUT);
852 
853  const char * telname = gravi_conf_get_telname (0, header);
854  CPLCHECK ("Cannot get telescope name");
855  hdrl_parameter * strehl_params;
856 
857  /* Hardcoded Mirror diameters */
858  if (telname[0] == 'U') {
859  strehl_params =
860  hdrl_strehl_parameter_create (1.8e-6, 8.0/2.0, 1.126/2.0,
861  pscale*1e-3, pscale*1e-3,
862  0.8, 0.8, 1.0);
863  } else if (telname[0] == 'A') {
864  strehl_params =
865  hdrl_strehl_parameter_create (1.8e-6, 1.8/2.0, 0.138/2.0,
866  pscale*1e-3, pscale*1e-3,
867  0.8, 0.8, 1.0);
868  } else {
869  cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
870  "Cannot recognise the telescope");
871  return CPL_ERROR_ILLEGAL_INPUT;
872  }
873 
874  hdrl_strehl_result strehl;
875  cpl_image * sub_image = cpl_image_extract (img, x-50, y-50, x+50, y+50);
876  hdrl_image * sub_hdrl = hdrl_image_create (sub_image, NULL);
877 
878  strehl = hdrl_strehl_compute (sub_hdrl, strehl_params);
879  *SR = (double) strehl.strehl_value.data;
880  FREE (hdrl_image_delete, sub_hdrl);
881  FREE (cpl_image_delete, sub_image);
882  hdrl_parameter_delete(strehl_params);
883 
884  gravi_msg_function_exit(0);
885  return CPL_ERROR_NONE;
886 }
887 
888 /*----------------------------------------------------------------------------*/
889 /*
890  * @brief measure maximum of the source at the given location
891  */
892 /*----------------------------------------------------------------------------*/
893 
894 cpl_error_code gravi_acq_measure_max(cpl_image * img, double x, double y, double size, double * img_max)
895 {
896  gravi_msg_function_start(0);
897  cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
898  cpl_ensure_code (img_max, CPL_ERROR_NULL_INPUT);
899 
900  cpl_image * small_img = cpl_image_extract(img, x-size, y-size, x+size, y+size);
901  *img_max = cpl_image_get_max(small_img);
902  cpl_image_delete(small_img);
903 
904  gravi_msg_function_exit(0);
905  return CPL_ERROR_NONE;
906 }
907 
908 
909 /*----------------------------------------------------------------------------*/
923 /*----------------------------------------------------------------------------*/
924 
925 cpl_error_code gravi_acq_fit_gaussian (cpl_image * img, double *x, double *y,
926  double *ex, double *ey, cpl_size size)
927 {
928  gravi_msg_function_start(0);
929  cpl_ensure_code (img, CPL_ERROR_NULL_INPUT);
930  cpl_ensure_code (x, CPL_ERROR_NULL_INPUT);
931  cpl_ensure_code (y, CPL_ERROR_NULL_INPUT);
932  cpl_ensure_code (ex, CPL_ERROR_NULL_INPUT);
933  cpl_ensure_code (ey, CPL_ERROR_NULL_INPUT);
934 
935  /* Fill first guess */
936  cpl_array * parameters = cpl_array_new (7, CPL_TYPE_DOUBLE);
937  cpl_array_fill_window_invalid (parameters, 0, 7);
938  cpl_array_set (parameters, 2, 0);
939  cpl_array_set (parameters, 3, *x);
940  cpl_array_set (parameters, 4, *y);
941  cpl_array_set (parameters, 5, 3);
942  cpl_array_set (parameters, 6, 3);
943  CPLCHECK_MSG ("Error creating parameter table");
944 
945  double med = cpl_image_get_median_window (img,
946  (cpl_size)(*x)-size, (cpl_size)(*y)-size,
947  (cpl_size)(*x)+size, (cpl_size)(*y)+size);
948  cpl_array_set (parameters, 0, med);
949  CPLCHECK_MSG ("Error getting median");
950 
951  /* Fit Gaussian */
952  double rms=0.;
953  cpl_fit_image_gaussian (img, NULL, (cpl_size)(*x), (cpl_size)(*y),
954  size, size, parameters,
955  NULL, NULL, &rms, NULL, NULL,
956  NULL, NULL, NULL, NULL);
957  CPLCHECK_MSG ("Error fitting Gaussian");
958 
959  /* Set back */
960  *x = cpl_array_get (parameters,3,NULL);
961  *y = cpl_array_get (parameters,4,NULL);
962  CPLCHECK_MSG ("Error reading fit result");
963 
964  /* Check errors */
965  /* reject result if either: */
966  /* * peak is below 3*rms */
967  /* * sigma_x > 5 pixels */
968  /* * sigma_y > 5 pixels */
969  double A = cpl_array_get (parameters, 1, NULL);
970  double rho = cpl_array_get (parameters, 2, NULL);
971  double sx = cpl_array_get (parameters, 5, NULL);
972  double sy = cpl_array_get (parameters, 6, NULL);
973 
974  if ( A < 0. ) {
975  // detection is just not significant
976  cpl_msg_info (cpl_func, "rejecting fit: x=%g, y=%g, SNR=%g, sx=%g, sy=%g",
977  *x, *y, A/(rms * 2.*M_PI*sx*sy*sqrt(1-rho*rho)), sx, sy);
978  *x = 0.;
979  *y = 0.;
980  *ex = -1.;
981  *ey = -1.;
982  } else {
983  // cf. Condon 1996, PASP 109:166
984  double cst = 2. * sqrt(2. * M_PI * (1. - rho*rho) * sx * sy) * rms / A;
985  *ex=cst*sx;
986  *ey=cst*sy;
987  }
988  CPLCHECK_MSG("Error checking significance of fit result");
989 
990  /* Fill image with zero at the detected position */
991  //if (*x > 0. && *y > 0.) {
992  // cpl_image_set (img, (cpl_size)(*x), (cpl_size)(*y), 0.0);
993  //}
994  //CPLCHECK_MSG("Error setting peak to zero");
995 
996  /* Delete */
997  FREE (cpl_array_delete, parameters);
998 
999  gravi_msg_function_exit(0);
1000  return CPL_ERROR_NONE;
1001 }
1002 
1003 /*----------------------------------------------------------------------------*/
1018 /*----------------------------------------------------------------------------*/
1019 
1020 cpl_error_code gravi_acqcam_pupil (cpl_image * mean_img,
1021  cpl_imagelist * acqcam_imglist,
1022  cpl_propertylist * header,
1023  cpl_table * acqcam_table,
1024  cpl_propertylist * o_header,
1025  gravi_data *static_param_data)
1026 {
1027  gravi_msg_function_start(1);
1028  cpl_ensure_code (mean_img, CPL_ERROR_NULL_INPUT);
1029  cpl_ensure_code (acqcam_imglist, CPL_ERROR_NULL_INPUT);
1030  cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1031  cpl_ensure_code (acqcam_table, CPL_ERROR_NULL_INPUT);
1032  cpl_ensure_code (o_header, CPL_ERROR_NULL_INPUT);
1033 
1034  /* Number of row */
1035  cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
1036 
1037  /* Compute separation */
1038  double sobj_x = gravi_pfits_get_sobj_x (header);
1039  double sobj_y = gravi_pfits_get_sobj_y (header);
1040  double sobj_sep = sqrt(sobj_x*sobj_x + sobj_y*sobj_y);
1041 
1042  /* Pupil positions (or we use array of 3) */
1043  gravi_table_new_column (acqcam_table, "PUPIL_NSPOT", NULL, CPL_TYPE_INT);
1044  gravi_table_new_column (acqcam_table, "PUPIL_X", "pix", CPL_TYPE_DOUBLE);
1045  gravi_table_new_column (acqcam_table, "PUPIL_Y", "pix", CPL_TYPE_DOUBLE);
1046  gravi_table_new_column (acqcam_table, "PUPIL_Z", "pix", CPL_TYPE_DOUBLE);
1047  gravi_table_new_column (acqcam_table, "PUPIL_R", "deg", CPL_TYPE_DOUBLE);
1048  gravi_table_new_column (acqcam_table, "PUPIL_U", "m", CPL_TYPE_DOUBLE);
1049  gravi_table_new_column (acqcam_table, "PUPIL_V", "m", CPL_TYPE_DOUBLE);
1050  gravi_table_new_column (acqcam_table, "PUPIL_W", "m", CPL_TYPE_DOUBLE);
1051  gravi_table_new_column (acqcam_table, "OPD_PUPIL", "m", CPL_TYPE_DOUBLE);
1052 
1053  int nspot = 0, ntel = 4;
1054 
1055  /* Loop on tel */
1056  for (int tel = 0; tel < ntel; tel++) {
1057  cpl_msg_info (cpl_func, "Compute pupil position for beam %i", tel+1);
1058 
1059  /* Get the conversion angle xy to uv in [rad] */
1060  double fangle = gravi_pfits_get_fangle_acqcam (header, tel);
1061  double cfangle = cos(fangle * CPL_MATH_RAD_DEG);
1062  double sfangle = sin(fangle * CPL_MATH_RAD_DEG);
1063  CPLCHECK_MSG ("Cannot read ESO INS DROTOFF#");
1064 
1065  /* Get the orientation of star */
1066  double drotoff = gravi_pfits_get_drotoff (header, tel);
1067  double cdrotoff = cos(drotoff * CPL_MATH_RAD_DEG);
1068  double sdrotoff = sin(drotoff * CPL_MATH_RAD_DEG);
1069  CPLCHECK_MSG ("Cannot read ESO INS DROTOFF#");
1070 
1071  /* Allocate memory */
1072  cpl_vector * a_start = cpl_vector_new (GRAVI_SPOT_NA);
1073  cpl_vector_fill (a_start, 0.0);
1074 
1075  /* Read the sub-apperture reference positions
1076  * Converted to accound for sub-windowing
1077  * In vector convention (start at 0,0) */
1078  gravi_acqcam_get_pup_ref (header, tel, a_start);
1079  gravi_acqcam_get_diode_ref (header, tel, a_start);
1080  CPLCHECK_MSG ("Cannot read ACQ data for pupil in header");
1081 
1082  cpl_vector * a_final = cpl_vector_duplicate (a_start);
1083  CPLCHECK_MSG ("Cannot prepare parameters");
1084 
1085  /* Fit pupil sensor spots in this image, with various
1086  * starting points to converge to the global minimum */
1087  gravi_acqcam_fit_spot (mean_img, 30, a_final, 1, &nspot);
1088  CPLCHECK_MSG ("Cannot fit rotation and center");
1089 
1090  /* Offsets in pixels */
1091  double xpos = cpl_vector_get (a_final,GRAVI_SPOT_SUB+0) -
1092  cpl_vector_get (a_start,GRAVI_SPOT_SUB+0);
1093  double ypos = cpl_vector_get (a_final,GRAVI_SPOT_SUB+4) -
1094  cpl_vector_get (a_start,GRAVI_SPOT_SUB+4);
1095 
1096  /* Scale and rotation in deg (rectangle, thus 180deg symetrie) */
1097  double scale = cpl_vector_get (a_final,GRAVI_SPOT_SCALE);
1098  double angle = cpl_vector_get (a_final,GRAVI_SPOT_ANGLE);
1099  if (angle < 0) angle += 180;
1100  if (angle > 180) angle -= 180;
1101  cpl_vector_set (a_final,GRAVI_SPOT_ANGLE,angle);
1102 
1103  /* In UV [m] */
1104  double upos = (cfangle * xpos - sfangle * ypos) / scale;
1105  double vpos = (sfangle * xpos + cfangle * ypos) / scale;
1106 
1107  /* Add best position as a cross in image */
1108  gravi_acqcam_spot_imprint (mean_img, a_final);
1109 
1110  /* Add QC parameters */
1111  char qc_name[100];
1112 
1113  sprintf (qc_name, "ESO QC ACQ FIELD%i NORTH_ANGLE", tel+1);
1114  cpl_msg_info (cpl_func, "%s = %f", qc_name, fangle);
1115  cpl_propertylist_update_double (o_header, qc_name, fangle);
1116  cpl_propertylist_set_comment (o_header, qc_name, "[deg] y->x, predicted North direction on ACQ");
1117 
1118  sprintf (qc_name, "ESO QC ACQ PUP%i NSPOT", tel+1);
1119  cpl_msg_info (cpl_func, "%s = %i", qc_name, nspot);
1120  cpl_propertylist_update_int (o_header, qc_name, nspot);
1121  cpl_propertylist_set_comment (o_header, qc_name, "nb. of pupil spot in ACQ");
1122 
1123  sprintf (qc_name, "ESO QC ACQ PUP%i ANGLE", tel+1);
1124  cpl_msg_info (cpl_func, "%s = %f", qc_name, angle);
1125  cpl_propertylist_update_double (o_header, qc_name, angle);
1126  cpl_propertylist_set_comment (o_header, qc_name, "[deg] y->x, diode angle on ACQ");
1127 
1128  sprintf (qc_name, "ESO QC ACQ PUP%i SCALE", tel+1);
1129  cpl_msg_info (cpl_func, "%s = %f", qc_name, scale);
1130  cpl_propertylist_update_double (o_header, qc_name, scale);
1131  cpl_propertylist_set_comment (o_header, qc_name, "[pix/m] diode scale on ACQ");
1132 
1133  sprintf (qc_name, "ESO QC ACQ PUP%i XPOS", tel+1);
1134  cpl_msg_info (cpl_func, "%s = %f", qc_name, xpos);
1135  cpl_propertylist_update_double (o_header, qc_name, xpos);
1136  cpl_propertylist_set_comment (o_header, qc_name, "[pix] pupil x-shift in ACQ");
1137 
1138  sprintf (qc_name, "ESO QC ACQ PUP%i YPOS", tel+1);
1139  cpl_msg_info (cpl_func, "%s = %f", qc_name, ypos);
1140  cpl_propertylist_update_double (o_header, qc_name, ypos);
1141  cpl_propertylist_set_comment (o_header, qc_name, "[pix] pupil y-shift in ACQ");
1142 
1143  sprintf (qc_name, "ESO QC ACQ PUP%i UPOS", tel+1);
1144  cpl_msg_info (cpl_func, "%s = %f", qc_name, upos);
1145  cpl_propertylist_update_double (o_header, qc_name, upos);
1146  cpl_propertylist_set_comment (o_header, qc_name, "[m] pupil u-shift in ACQ");
1147 
1148  sprintf (qc_name, "ESO QC ACQ PUP%i VPOS", tel+1);
1149  cpl_msg_info (cpl_func, "%s = %f", qc_name, vpos);
1150  cpl_propertylist_update_double (o_header, qc_name, vpos);
1151  cpl_propertylist_set_comment (o_header, qc_name, "[m] pupil v-shift in ACQ");
1152 
1153  /* Loop on all images */
1154  for (cpl_size row = 0; row < nrow; row++) {
1155  if (row %10 == 0 || row == (nrow-1))
1156  cpl_msg_info_overwritable (cpl_func, "Fit image %lld over %lld", row+1, nrow);
1157 
1158  /* Get data and first guess */
1159  cpl_image * img = cpl_imagelist_get (acqcam_imglist, row);
1160  cpl_vector * a_row = cpl_vector_duplicate (a_final);
1161 
1162  /* Fit pupil sensor spots in this image. */
1163  gravi_acqcam_fit_spot (img, 1, a_row, 0, &nspot);
1164  CPLCHECK_MSG ("Cannot fit sub-appertures of image");
1165 
1166  /* If spot detected */
1167  if (nspot < 4) {
1168  cpl_table_set (acqcam_table, "PUPIL_NSPOT", row*ntel+tel, 0);
1169  }
1170  else {
1171  /* Add best position as a cross in image */
1172  gravi_acqcam_spot_imprint (img, a_row);
1173 
1174  /* Remove reference */
1175  cpl_vector_subtract (a_row, a_start);
1176 
1177  /* Get rotation [deg] and shift [pix]*/
1178  double r_shift = cpl_vector_get (a_row, GRAVI_SPOT_ANGLE);
1179  double x_shift = cpl_vector_get (a_row, GRAVI_SPOT_SUB+0);
1180  double y_shift = cpl_vector_get (a_row, GRAVI_SPOT_SUB+4);
1181  double z_shift = -0.5 * ( cpl_vector_get (a_row, GRAVI_SPOT_SUB+2) +
1182  cpl_vector_get (a_row, GRAVI_SPOT_SUB+5));
1183 
1184  /* In UV [m] */
1185  double u_shift = (cfangle * x_shift - sfangle * y_shift) / scale;
1186  double v_shift = (sfangle * x_shift + cfangle * y_shift) / scale;
1187  double w_shift = gravi_acqcam_z2meter (z_shift, static_param_data);
1188 
1189  cpl_table_set (acqcam_table, "PUPIL_NSPOT", row*ntel+tel, nspot);
1190  cpl_table_set (acqcam_table, "PUPIL_X", row*ntel+tel, x_shift);
1191  cpl_table_set (acqcam_table, "PUPIL_Y", row*ntel+tel, y_shift);
1192  cpl_table_set (acqcam_table, "PUPIL_Z", row*ntel+tel, z_shift);
1193  cpl_table_set (acqcam_table, "PUPIL_R", row*ntel+tel, r_shift);
1194  cpl_table_set (acqcam_table, "PUPIL_U", row*ntel+tel, u_shift);
1195  cpl_table_set (acqcam_table, "PUPIL_V", row*ntel+tel, v_shift);
1196  cpl_table_set (acqcam_table, "PUPIL_W", row*ntel+tel, w_shift);
1197 
1198  /* Compute the OPD_PUPIL */
1199  double opd_pupil = sobj_sep * GRAVI_MATH_RAD_MAS / scale *
1200  (x_shift * cdrotoff + y_shift * sdrotoff);
1201  cpl_table_set (acqcam_table, "OPD_PUPIL", row*ntel+tel, opd_pupil);
1202  }
1203 
1204  FREE (cpl_vector_delete, a_row);
1205  }
1206 
1207  FREE (cpl_vector_delete, a_start);
1208  FREE (cpl_vector_delete, a_final);
1209  } /* End loop on tel */
1210 
1211  gravi_msg_function_exit(1);
1212  return CPL_ERROR_NONE;
1213 }
1214 
1215 /*----------------------------------------------------------------------------*/
1232 /*----------------------------------------------------------------------------*/
1233 
1234 cpl_error_code gravi_acqcam_field (cpl_image * mean_img,
1235  cpl_imagelist * acqcam_imglist,
1236  cpl_propertylist * header,
1237  cpl_table * acqcam_table,
1238  cpl_propertylist * o_header,
1239  gravi_data *static_param_data)
1240 {
1241  gravi_msg_function_start(1);
1242  cpl_ensure_code (mean_img, CPL_ERROR_NULL_INPUT);
1243  cpl_ensure_code (acqcam_imglist, CPL_ERROR_NULL_INPUT);
1244  cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
1245  cpl_ensure_code (acqcam_table, CPL_ERROR_NULL_INPUT);
1246  cpl_ensure_code (o_header, CPL_ERROR_NULL_INPUT);
1247 
1248  char qc_name[100];
1249  int ntel = 4;
1250 
1251  /* Number of row */
1252  cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
1253 
1254  /* Create columns */
1255  gravi_table_new_column (acqcam_table, "FIELD_SC_X", "pix", CPL_TYPE_DOUBLE);
1256  gravi_table_new_column (acqcam_table, "FIELD_SC_Y", "pix", CPL_TYPE_DOUBLE);
1257  gravi_table_new_column (acqcam_table, "FIELD_SC_XERR", "pix", CPL_TYPE_DOUBLE);
1258  gravi_table_new_column (acqcam_table, "FIELD_SC_YERR", "pix", CPL_TYPE_DOUBLE);
1259  gravi_table_new_column (acqcam_table, "FIELD_FT_X", "pix", CPL_TYPE_DOUBLE);
1260  gravi_table_new_column (acqcam_table, "FIELD_FT_Y", "pix", CPL_TYPE_DOUBLE);
1261  gravi_table_new_column (acqcam_table, "FIELD_FT_XERR", "pix", CPL_TYPE_DOUBLE);
1262  gravi_table_new_column (acqcam_table, "FIELD_FT_YERR", "pix", CPL_TYPE_DOUBLE);
1263  gravi_table_new_column (acqcam_table, "FIELD_SCALE", "mas/pix", CPL_TYPE_DOUBLE);
1264  gravi_table_new_column (acqcam_table, "FIELD_SCALEERR", "mas/pix", CPL_TYPE_DOUBLE);
1265  gravi_table_new_column (acqcam_table, "FIELD_FIBER_DX", "pix", CPL_TYPE_DOUBLE);
1266  gravi_table_new_column (acqcam_table, "FIELD_FIBER_DY", "pix", CPL_TYPE_DOUBLE);
1267  gravi_table_new_column (acqcam_table, "FIELD_FIBER_DXERR", "pix", CPL_TYPE_DOUBLE);
1268  gravi_table_new_column (acqcam_table, "FIELD_FIBER_DYERR", "pix", CPL_TYPE_DOUBLE);
1269  gravi_table_new_column (acqcam_table, "FIELD_STREHL", "ratio", CPL_TYPE_DOUBLE);
1270 
1271  /* ----------------- START EKW 26/11/2018 read constant parameter from calibration file */
1272  /* Position of roof center on full frame */
1273  //double roof_x[] = {274.4, 787.1, 1236.1, 1673.4};
1274  //double roof_y[] = {242.3, 247.7, 225.8, 235.6};
1275  double *roof_x = NULL;
1276  double *roof_y = NULL;
1277 
1278  /* Position of single-field spot on full frame */
1279  //double spot_x[] = {289. , 798.2, 1245.5, 1696.};
1280  //double spot_y[] = {186.5, 187.5, 172.5, 178.};
1281  double *spot_x = NULL;
1282  double *spot_y = NULL;
1283 
1284  /* Default position angle of roof */
1285  //double roof_pos[] = {38.49, 38.54, 38.76, 39.80};
1286  double *roof_pos = NULL;
1287 
1288  cpl_table * roof_pos_table = gravi_data_get_table (static_param_data, "ROOFPOS");
1289  CPLCHECK_MSG ("STATIC_CALIB not available in the SOF. It is mandatory for acqcam reduction.");
1290 
1291  if ( cpl_table_has_column(roof_pos_table , "roof_x") ) {
1292  roof_x= cpl_table_get_data_double (roof_pos_table, "roof_x");
1293  cpl_msg_info(cpl_func,"roof_x [0] : %e", roof_x[0] );
1294  cpl_msg_info(cpl_func,"roof_x [1] : %e", roof_x[1] );
1295  cpl_msg_info(cpl_func,"roof_x [2] : %e", roof_x[2] );
1296  cpl_msg_info(cpl_func,"roof_x [3] : %e", roof_x[3] );
1297  }
1298  else {
1299  cpl_msg_warning(cpl_func,"Cannot get the default values for roof_x ");
1300  }
1301 
1302  if ( cpl_table_has_column(roof_pos_table , "roof_y") ) {
1303  roof_y= cpl_table_get_data_double (roof_pos_table, "roof_y");
1304  cpl_msg_info(cpl_func,"roof_y [0] : %e", roof_y[0] );
1305  cpl_msg_info(cpl_func,"roof_y [1] : %e", roof_y[1] );
1306  cpl_msg_info(cpl_func,"roof_y [2] : %e", roof_y[2] );
1307  cpl_msg_info(cpl_func,"roof_y [3] : %e", roof_y[3] );
1308  }
1309  else {
1310  cpl_msg_warning(cpl_func,"Cannot get the default values for roof_y ");
1311  }
1312 
1313  if ( cpl_table_has_column(roof_pos_table , "spot_x") ) {
1314  spot_x= cpl_table_get_data_double (roof_pos_table, "spot_x");
1315  cpl_msg_info(cpl_func,"spot_x [0] : %e", spot_x[0] );
1316  cpl_msg_info(cpl_func,"spot_x [1] : %e", spot_x[1] );
1317  cpl_msg_info(cpl_func,"spot_x [2] : %e", spot_x[2] );
1318  cpl_msg_info(cpl_func,"spot_x [3] : %e", spot_x[3] );
1319  }
1320  else {
1321  cpl_msg_warning(cpl_func,"Cannot get the default values for spot_x ");
1322  }
1323 
1324  if ( cpl_table_has_column(roof_pos_table , "spot_y") ) {
1325  spot_y= cpl_table_get_data_double (roof_pos_table, "spot_y");
1326  cpl_msg_info(cpl_func,"spot_y [0] : %e", spot_y[0] );
1327  cpl_msg_info(cpl_func,"spot_y [1] : %e", spot_y[1] );
1328  cpl_msg_info(cpl_func,"spot_y [2] : %e", spot_y[2] );
1329  cpl_msg_info(cpl_func,"spot_y [3] : %e", spot_y[3] );
1330  }
1331  else {
1332  cpl_msg_warning(cpl_func,"Cannot get the default values for spot_y ");
1333  }
1334 
1335  if ( cpl_table_has_column(roof_pos_table , "roof_pos") ) {
1336  roof_pos= cpl_table_get_data_double (roof_pos_table, "roof_pos");
1337  cpl_msg_info(cpl_func,"roof_pos [0] : %e", roof_pos[0] );
1338  cpl_msg_info(cpl_func,"roof_pos [1] : %e", roof_pos[1] );
1339  cpl_msg_info(cpl_func,"roof_pos [2] : %e", roof_pos[2] );
1340  cpl_msg_info(cpl_func,"roof_pos [3] : %e", roof_pos[3] );
1341  }
1342  else {
1343  cpl_msg_warning(cpl_func,"Cannot get the default values for roof_pos ");
1344  }
1345  /* ------------------ END EKW 26/11/2018 read constant parameter from calibration file */
1346 
1347 
1348 
1349  /* Static acqcam field calibration error */
1350  /* calibration from Frank, March 2017, using Marcel
1351  double calib_dx[] = {-0.3483, -1.0251, -0.5432, -0.2024} ;
1352  double calib_dy[] = { 0.3089, -0.5221, -0.2686, -0.3843} ; */
1353  /* calibration from Julien, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-10
1354  double calib_dx[] = {-0.61, -1.69, -0.97, 0.00} ;
1355  double calib_dy[] = { 0.53, -1.07, -0.49, -0.47} ; */
1356  /* calibration from Oli, using GJ65 observation with Sylvestre tracking on 2018-01-03, email 2018-01-19/
1357  double calib_dx[] = {-0.63, -1.67, -0.97, 0.02} ;
1358  double calib_dy[] = { 0.55, -1.07, -0.49, -0.45} ; */
1359  /* FE: removing obsolete calib_dx, because now handled by ICS, or by
1360  correcting fitsheader retrospectively */
1361  double calib_dx[] = {0,0,0,0} ;
1362  double calib_dy[] = {0,0,0,0} ;
1363 
1364  /* If sub-windowing, we read the sub-window size */
1365  cpl_size nsx = 512;
1366  /* EKW 10/01/2019 cpl_size nsy = 512; */
1367  if ( cpl_propertylist_has (header, "ESO DET1 FRAMES NX") ) {
1368  nsx = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NX");
1369  /* EKW 10/01/2019 nsy = cpl_propertylist_get_int (header, "ESO DET1 FRAMES NY"); */
1370  }
1371 
1372  /* Compute separation */
1373  double sobj_x = gravi_pfits_get_sobj_x (header);
1374  double sobj_y = gravi_pfits_get_sobj_y (header);
1375 
1376  /* Accumulated mapping offsets */
1377  double sobj_offx=0., sobj_offy=0, sobj_drho=0., sobj_dth=0.;
1378  if (cpl_propertylist_has(header, "ESO INS SOBJ OFFX")) {
1379  /* mapping / mosaicing blind offsets from otiginal acquisition
1380  * to current position */
1381  sobj_offx = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFX");
1382  sobj_offy = cpl_propertylist_get_double (header, "ESO INS SOBJ OFFY");
1383 
1384  /* distance from acquired SC position to current SC position in mas */
1385  sobj_drho = sqrt(sobj_offx*sobj_offx+sobj_offy*sobj_offy);
1386 
1387  /* position angle on sky from acquired SC position to */
1388  /* current SC position, neglecting anamorphism variations */
1389  sobj_dth = atan2(sobj_offx, sobj_offy)*180./M_PI;
1390  }
1391  CPLCHECK_MSG ("Cannot get separation");
1392 
1393  /* Recover position of originally acquired star, before any blind offset */
1394  double dx_in = sobj_x - sobj_offx;
1395  double dy_in = sobj_y - sobj_offy;
1396  double rho_in = sqrt(dx_in*dx_in + dy_in*dy_in);
1397 
1398  /* Force zero separation for SINGLE mode */
1399  if ( gravi_pfits_get_mode (header) == MODE_SINGLE) rho_in = 0.;
1400  CPLCHECK_MSG ("Error reading header information");
1401 
1402  /* Loop on tel */
1403  for (int tel = 0; tel < ntel; tel++) {
1404  char name[90];
1405  double rp=roof_pos[tel]; // default value of roof position angle
1406  cpl_size sx=tel*512+1, sy=1;
1407  cpl_msg_info (cpl_func, "Compute field position for beam %i", tel+1);
1408 
1409  /* Read roof position angle */
1410  sprintf (name, "ESO INS DROTOFF%d", tel + 1);
1411  if ( cpl_propertylist_has (header, name) ) {
1412  rp = cpl_propertylist_get_double(header, name);
1413  CPLCHECK ("Cannot get rotation");
1414  }
1415 
1416  /* Approx. position angle of the binary, left from top */
1417  double approx_PA = 270.-rp;
1418 
1419  /* Get the telescope name and ID */
1420  const char * telname = gravi_conf_get_telname (tel, header);
1421 
1422  /* Check telescope name */
1423  if (!telname) cpl_msg_error (cpl_func, "Cannot read the telescope name");
1424  cpl_ensure_code (telname, CPL_ERROR_ILLEGAL_INPUT);
1425 
1426  /* Hardcoded approx. plate-scale in mas/pix. */
1427  double scale = 0.0;
1428  char telid[128];
1429  if (telname[0] == 'U') {
1430  /* scale = 18.; */
1431  scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO PLATE SCALE UT");
1432  cpl_msg_info (cpl_func,"PLATE SCALE UT is : %e \n",scale);
1433  } else if (telname[0] == 'A') {
1434  sprintf(telid, "ESO PLATE SCALE AT%d", tel+1);
1435  scale = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), telid);
1436  cpl_msg_info (cpl_func,"PLATE SCALE AT%d is : %e \n",tel+1, scale);
1437  /*
1438  if (tel == 0) {
1439  scale = 76.8;
1440  } else if (tel == 1) {
1441  scale = 78.0;
1442  } else if (tel == 2) {
1443  scale = 77.0;
1444  } else if (tel == 3) {
1445  scale = 84.6;
1446  }
1447  */
1448  }
1449 
1450  /* Position of the fibres */
1451  double fiber_xft=0.;
1452  double fiber_yft=0.;
1453  double fiber_xsc=0.;
1454  double fiber_ysc=0.;
1455  sprintf (name, "ESO ACQ FIBER FT%dX", tel + 1);
1456  if (cpl_propertylist_has(header, name)) {
1457  fiber_xft=cpl_propertylist_get_double(header, name);
1458  sprintf (name, "ESO ACQ FIBER FT%dY", tel + 1);
1459  fiber_yft=cpl_propertylist_get_double(header, name);
1460  sprintf (name, "ESO ACQ FIBER SC%dX", tel + 1);
1461  fiber_xsc=cpl_propertylist_get_double(header, name);
1462  sprintf (name, "ESO ACQ FIBER SC%dY", tel + 1);
1463  fiber_ysc=cpl_propertylist_get_double(header, name);
1464  }
1465  double fiber_ft_sc_x=fiber_xsc-fiber_xft;
1466  double fiber_ft_sc_y=fiber_ysc-fiber_yft;
1467 
1468  /* Get the North position angle on the camera */
1469  double fangle = gravi_pfits_get_fangle_acqcam (header, tel);
1470  CPLCHECK ("Cannot get rotation");
1471 
1472  /* Mapping/mosaicing offset on acq cam axes, in mas, */
1473  /* neglecting amnamorphism variations */
1474  double sobj_offx_cam = sobj_drho * sin((fangle+sobj_dth)/180.*M_PI);
1475  double sobj_offy_cam = sobj_drho * cos((fangle+sobj_dth)/180.*M_PI);
1476 
1477  /* If sub-windowing, we read the sub-window start for field */
1478  if ( nsx != 512 ) {
1479  sprintf (name, "ESO DET1 FRAM%d STRX", tel + 1);
1480  sx = cpl_propertylist_get_int (header, name);
1481 
1482  sprintf (name, "ESO DET1 FRAM%d STRY", tel + 1);
1483  sy = cpl_propertylist_get_int (header, name);
1484 
1485  CPLCHECK_MSG ("Cannot get sub-windowing parameters");
1486  }
1487 
1488  cpl_msg_debug (cpl_func,"sub-window field %i sx= %lld sy = %lld", tel, sx, sy);
1489 
1490  /*--------------------------------------------------*/
1491  /* Guess position of SC and FT targes in mean image */
1492  /*--------------------------------------------------*/
1493 
1494  double xFT, yFT, xSC, ySC;
1495 
1496  if (rho_in == 0.) {
1497  /* TODO: close dual-field */
1498  /* Single-field case */
1499  /* Simply shift the best spot from full frame to cut-out */
1500  xFT = spot_x[tel] - sx + nsx*tel + 1;
1501  yFT = spot_y[tel] - sy + 1;
1502  xSC = xFT;
1503  ySC = yFT;
1504  } else {
1505  /* Pixel position of roof center on cut-out frame */
1506  /* Shift from full frame to cut-out */
1507  double cutout_roof_x = roof_x[tel] - sx + nsx*tel + 1;
1508  double cutout_roof_y = roof_y[tel] - sy + 1;
1509 
1510  cpl_msg_info (cpl_func, " ROOF_X = %f", cutout_roof_x);
1511  cpl_msg_info (cpl_func, " ROOF_Y = %f", cutout_roof_y);
1512 
1513  /* Approx pixel offset from SC to FT, divided by 2 */
1514  double approx_dx=0.5*rho_in*sin(approx_PA*M_PI/180.)/scale;
1515  double approx_dy=0.5*rho_in*cos(approx_PA*M_PI/180.)/scale;
1516 
1517  /* Expected position of the two stars */
1518  xFT = cutout_roof_x - approx_dx ;
1519  yFT = cutout_roof_y - approx_dy ;
1520  xSC = cutout_roof_x + approx_dx ;
1521  ySC = cutout_roof_y + approx_dy ;
1522  }
1523 
1524  cpl_msg_info (cpl_func, "guess SC_X = %f", xSC);
1525  cpl_msg_info (cpl_func, "guess SC_Y = %f", ySC);
1526  cpl_msg_info (cpl_func, "guess FT_X = %f", xFT);
1527  cpl_msg_info (cpl_func, "guess FT_Y = %f", yFT);
1528 
1529  /*------------------------------------------*/
1530  /* Measure SC target position in mean image */
1531  /*------------------------------------------*/
1532 
1533  /* Box size */
1534  /* Optimal size has been determined empirically. */
1535  /* Too small value (15) will sometimes miss even a bright target. */
1536  /* Too large value (25) will often pick a wrong star (or backround noise) */
1537  cpl_size size = 20;
1538 
1539  double xSCguess=xSC, ySCguess=ySC, xFTguess=xFT, yFTguess=yFT;
1540  double ex, ey;
1541  double qc_val=0.;
1542 
1543  gravi_acq_fit_gaussian (mean_img, &xSC, &ySC, &ex, &ey, size);
1544  CPLCHECK_MSG("Error fitting SC");
1545 
1546  /*--------------------------------------------------------*/
1547  /* Add QC parameters for SC target position in mean image */
1548  /*--------------------------------------------------------*/
1549 
1550  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X", tel+1);
1551  if (xSC==0.) {
1552  /* Fitting failed: put QC to 0., reset xSC to gues value */
1553  qc_val = 0.;
1554  xSC = xSCguess;
1555  } else {
1556  /* Fiting succeeded: shift into full frame */
1557  qc_val = xSC + sx - 1 - nsx*tel;
1558  }
1559  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1560  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1561  cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
1562 
1563  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y", tel+1);
1564  if (ySC==0.) {
1565  /* Fitting failed: put QC to 0., reset ySC to gues value */
1566  qc_val = 0.;
1567  ySC = ySCguess;
1568  } else {
1569  /* Fiting succeeded: shift into full frame */
1570  qc_val = ySC + sy -1;
1571  }
1572  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1573  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1574  cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
1575 
1576  /*------------------------------------------*/
1577  /* Measure FT target position in mean image */
1578  /*------------------------------------------*/
1579 
1580  if (rho_in != 0.) {
1581  /* Detec FT */
1582  gravi_acq_fit_gaussian (mean_img, &xFT, &yFT, &ex, &ey, size);
1583  CPLCHECK_MSG("Error fitting FT");
1584  } else {
1585  xFT=xSC;
1586  yFT=ySC;
1587  }
1588 
1589  /*--------------------------------------------------------*/
1590  /* Add QC parameters for SC target position in mean image */
1591  /*--------------------------------------------------------*/
1592 
1593  sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X", tel+1);
1594  if (xFT==0.) {
1595  /* Fitting failed: put QC to 0., reset xFT to gues value */
1596  qc_val = 0.;
1597  xFT = xFTguess;
1598  } else {
1599  /* Fiting succeeded: shift into full frame */
1600  qc_val = xFT + sx - 1 - nsx*tel;
1601  }
1602  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1603  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1604  cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
1605 
1606  sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y", tel+1);
1607  if (yFT==0.) {
1608  /* Fitting failed: put QC to 0., reset yFT to gues value */
1609  qc_val = 0.;
1610  yFT = yFTguess;
1611  } else {
1612  /* Fiting succeeded: shift into full frame */
1613  qc_val = yFT + sy -1;
1614  }
1615  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1616  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1617  cpl_propertylist_set_comment (o_header, qc_name, "[pixel] position in mean image");
1618 
1619  /*---------------------------------------*/
1620  /* If in dual field, measure plate scale */
1621  /*---------------------------------------*/
1622 
1623  if (rho_in != 0.) {
1624  sprintf (qc_name, "ESO QC ACQ FIELD%i SCALE", tel+1);
1625  double sep = sqrt((ySC-yFT)*(ySC-yFT)+(xSC-xFT)*(xSC-xFT));
1626  double pscale = sep ? rho_in/sep : 0.;
1627  qc_val=pscale;
1628  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1629  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1630  cpl_propertylist_set_comment (o_header, qc_name,
1631  "[mas/pixel] plate-scale in the "
1632  "FT-SC direction");
1633  }
1634 
1635  /*------------------------------------------------*/
1636  /* If in dual field, measure fibre position error */
1637  /*------------------------------------------------*/
1638 
1639  if (rho_in != 0.) {
1640  /* Error in SC fibre positioning */
1641  /* The three terms are */
1642  /* - offset from FT target as detected to original SC */
1643  /* target as detected; */
1644  /* - blind offset command from mapping template projected */
1645  /* on acqisition camera; */
1646  /* - offset from FT fiber to SC fiber. */
1647  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DX", tel+1);
1648  qc_val = 0;
1649  if (scale) {
1650  qc_val = (xSC-xFT) + sobj_offx_cam/scale - fiber_ft_sc_x - calib_dx[tel];
1651  }
1652  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1653  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1654  cpl_propertylist_set_comment (o_header, qc_name,
1655  "[pixel] dx from SC fiber to "
1656  "SC object");
1657 
1658  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_FIBER_DY", tel+1);
1659  qc_val = 0;
1660  if (scale) {
1661  qc_val = (ySC-yFT) + sobj_offy_cam/scale - fiber_ft_sc_y - calib_dy[tel];
1662  }
1663  cpl_msg_info (cpl_func, "%s = %f", qc_name, qc_val);
1664  cpl_propertylist_update_double (o_header, qc_name, qc_val);
1665  cpl_propertylist_set_comment (o_header, qc_name,
1666  "[pixel] dx from SC fiber to "
1667  "SC object");
1668  }
1669 
1670  /*-----------------------------*/
1671  /* Average Strehl computation */
1672  /*-----------------------------*/
1673 
1674  double max_on_average, strehl_on_average;
1675 
1676  /* Measure strehl and maximum on averaged image */
1677  gravi_acq_measure_strehl(mean_img, xFT, yFT, scale, &strehl_on_average, header);
1678  gravi_acq_measure_max(mean_img, xFT, yFT, 15, &max_on_average);
1679 
1680  /* Update Strehl QC */
1681  sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL", tel+1);
1682  cpl_msg_info (cpl_func, "%s = %f", qc_name, strehl_on_average);
1683  cpl_propertylist_update_double (o_header, qc_name, strehl_on_average);
1684  cpl_propertylist_set_comment (o_header, qc_name, "Average Strehl from stacked AcqCam images");
1685 
1686  /*----------------------------*/
1687  /* Now process frame by frame */
1688  /*----------------------------*/
1689 
1690  /* Loop on all images */
1691  for (cpl_size row = 0; row < nrow; row++) {
1692  if (row %10 == 0 || row == (nrow-1))
1693  cpl_msg_info_overwritable (cpl_func, "Fit image %lld over %lld", row+1, nrow);
1694 
1695  /* Get data */
1696  cpl_image * img = cpl_imagelist_get (acqcam_imglist, row);
1697  CPLCHECK_MSG("Error getting image");
1698 
1699  /*-------------------------------------------------*/
1700  /* SC target position computation of current frame */
1701  /*-------------------------------------------------*/
1702  double xsc = xSC, ysc = ySC, exsc=0., eysc=0.;
1703  gravi_acq_fit_gaussian (img, &xsc, &ysc, &exsc, &eysc, size);
1704  CPLCHECK_MSG("Error fitting SC");
1705  /* Shift back positions to full frame */
1706  if (xsc != 0.) xsc += sx - 1 - nsx*tel;
1707  if (ysc != 0.) ysc += sy - 1;
1708  cpl_table_set (acqcam_table, "FIELD_SC_X", row*ntel+tel, xsc);
1709  cpl_table_set (acqcam_table, "FIELD_SC_Y", row*ntel+tel, ysc);
1710  cpl_table_set (acqcam_table, "FIELD_SC_XERR", row*ntel+tel, exsc);
1711  cpl_table_set (acqcam_table, "FIELD_SC_YERR", row*ntel+tel, eysc);
1712  CPLCHECK_MSG("Error setting SC columns");
1713 
1714  /*-------------------------------------------------*/
1715  /* FT target position computation of current frame */
1716  /*-------------------------------------------------*/
1717  double xft = xFT, yft = yFT, exft=0., eyft=0.;
1718  if (rho_in != 0.) {
1719  gravi_acq_fit_gaussian (img, &xft, &yft, &exft, &eyft, size);
1720  CPLCHECK_MSG("Error fitting FT");
1721  /* Shift back positions to full frame */
1722  if (xft != 0.) xft += sx - 1 - nsx*tel;
1723  if (yft != 0.) yft += sy - 1;
1724  } else {
1725  xft=xsc;
1726  yft=ysc;
1727  exft=exsc;
1728  eyft=eysc;
1729  }
1730  cpl_table_set (acqcam_table, "FIELD_FT_X", row*ntel+tel, xft);
1731  cpl_table_set (acqcam_table, "FIELD_FT_Y", row*ntel+tel, yft);
1732  cpl_table_set (acqcam_table, "FIELD_FT_XERR", row*ntel+tel, exft);
1733  cpl_table_set (acqcam_table, "FIELD_FT_YERR", row*ntel+tel, eyft);
1734  CPLCHECK_MSG("Error setting FT column");
1735 
1736  /*------------------------------------------*/
1737  /* Plate scale computation of current frame */
1738  /*------------------------------------------*/
1739  double ft_sc_x = xsc - xft;
1740  double ft_sc_y = ysc - yft;
1741  double eft_sc_x = sqrt(exsc*exsc+exft*exft);
1742  double eft_sc_y = sqrt(eysc*eysc+eyft*eyft);
1743  double sep = sqrt(ft_sc_x*ft_sc_x+ft_sc_y*ft_sc_y);
1744  double pscale = sep ? rho_in/sep : 0.;
1745  double escale = 0.;
1746  if (sep) {
1747  escale = rho_in/(sep*sep*sep)*
1748  (ft_sc_x*eft_sc_x+ft_sc_y*eft_sc_y) ;
1749  }
1750  cpl_table_set (acqcam_table, "FIELD_SCALE", row*ntel+tel, pscale);
1751  cpl_table_set (acqcam_table, "FIELD_SCALEERR", row*ntel+tel, escale);
1752 
1753  /*---------------------------------------------*/
1754  /* Error in fibre positioning of current frame */
1755  /*---------------------------------------------*/
1756  /* The three terms are */
1757  /* - offset from FT target as detected to original SC */
1758  /* target as detected; */
1759  /* - blind offset command from mapping template projected */
1760  /* on acqisition camera; */
1761  /* - offset from FT fiber to SC fiber. */
1762  double corrx=0, corry=0., ecorrx=0., ecorry=0.;
1763  if (pscale) {
1764  corrx = ft_sc_x + sobj_offx_cam/pscale - fiber_ft_sc_x - calib_dx[tel];
1765  corry = ft_sc_y + sobj_offy_cam/pscale - fiber_ft_sc_y - calib_dy[tel];
1766  double tmp = escale/(pscale*pscale);
1767  tmp *= tmp;
1768  ecorrx = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offx_cam*sobj_offx_cam*tmp);
1769  ecorry = sqrt(eft_sc_x*eft_sc_x + eft_sc_y*eft_sc_y + sobj_offy_cam*sobj_offy_cam*tmp);
1770  }
1771  cpl_table_set (acqcam_table, "FIELD_FIBER_DX", row*ntel+tel, corrx);
1772  cpl_table_set (acqcam_table, "FIELD_FIBER_DY", row*ntel+tel, corry);
1773  cpl_table_set (acqcam_table, "FIELD_FIBER_DXERR", row*ntel+tel, ecorrx);
1774  cpl_table_set (acqcam_table, "FIELD_FIBER_DYERR", row*ntel+tel, ecorry);
1775 
1776  /*-------------------------------------*/
1777  /* Strehl computation of current frame */
1778  /*-------------------------------------*/
1779  double max_on_frame;
1780  gravi_acq_measure_max(img, xFT, yFT, 15, &max_on_frame);
1781  cpl_table_set (acqcam_table, "FIELD_STREHL", row*ntel+tel, strehl_on_average*(max_on_frame/max_on_average) );
1782 
1783  } /* End loop on images */
1784 
1785  /* Add some QC */
1786  double sc_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_SC_X", tel, ntel);
1787  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_X STD", tel+1);
1788  cpl_propertylist_update_double (o_header, qc_name, sc_std_x);
1789  cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
1790 
1791  double sc_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_SC_Y", tel, ntel);
1792  sprintf (qc_name, "ESO QC ACQ FIELD%i SC_Y STD", tel+1);
1793  cpl_propertylist_update_double (o_header, qc_name, sc_std_y);
1794  cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of SC");
1795 
1796  double ft_std_x = gravi_table_get_column_std (acqcam_table, "FIELD_FT_X", tel, ntel);
1797  sprintf (qc_name, "ESO QC ACQ FIELD%i FT_X STD", tel+1);
1798  cpl_propertylist_update_double (o_header, qc_name, ft_std_x);
1799  cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
1800 
1801  double ft_std_y = gravi_table_get_column_std (acqcam_table, "FIELD_FT_Y", tel, ntel);
1802  sprintf (qc_name, "ESO QC ACQ FIELD%i FT_Y STD", tel+1);
1803  cpl_propertylist_update_double (o_header, qc_name, ft_std_y);
1804  cpl_propertylist_set_comment (o_header, qc_name, "[pix] Std of field position of FT");
1805 
1806  double strehl_std = gravi_table_get_column_std (acqcam_table, "FIELD_STREHL", tel, ntel);
1807  sprintf (qc_name, "ESO QC ACQ FIELD%i STREHL STD", tel+1);
1808  cpl_propertylist_update_double (o_header, qc_name, strehl_std);
1809  cpl_propertylist_set_comment (o_header, qc_name, "Std of FT strehl");
1810 
1811  } /* End loop on tel */
1812 
1813  gravi_msg_function_exit(1);
1814  return CPL_ERROR_NONE;
1815 }
1816 
1817 /*----------------------------------------------------------------------------*/
1834 /*----------------------------------------------------------------------------*/
1835 
1836 cpl_error_code gravi_reduce_acqcam (gravi_data * output_data,
1837  gravi_data * input_data,
1838  gravi_data * static_param_data)
1839 {
1840  gravi_msg_function_start(1);
1841  cpl_ensure_code (output_data, CPL_ERROR_NULL_INPUT);
1842  cpl_ensure_code (input_data, CPL_ERROR_NULL_INPUT);
1843 
1844  /* Check if extension exist */
1845  if (!gravi_data_has_extension (input_data, GRAVI_IMAGING_DATA_ACQ_EXT)) {
1846  gravi_msg_warning (cpl_func, "Cannot reduce the ACQCAM, not data");
1847  return CPL_ERROR_NONE;
1848  }
1849 
1850  /* Get the data and header */
1851  cpl_propertylist * header, * o_header;
1852  header = gravi_data_get_header (input_data);
1853  o_header = gravi_data_get_header (output_data);
1854 
1855  cpl_imagelist * acqcam_imglist;
1856  acqcam_imglist = gravi_data_get_cube (input_data, GRAVI_IMAGING_DATA_ACQ_EXT);
1857  CPLCHECK_MSG ("Cannot get data or header");
1858 
1859  /* Build the table */
1860  cpl_size ntel = 4;
1861  cpl_size nrow = cpl_imagelist_get_size (acqcam_imglist);
1862 
1863  cpl_table * acqcam_table;
1864  acqcam_table = cpl_table_new (nrow * ntel);
1865 
1866  /* Time column shall contain the time from PCR.ACQ.START in [us] */
1867  cpl_table_new_column (acqcam_table, "TIME", CPL_TYPE_INT);
1868  cpl_table_set_column_unit (acqcam_table, "TIME", "us");
1869 
1870  /* Loop on DIT in cube to fill the TIME column
1871  * same value for all 4 beams*/
1872  for (cpl_size row = 0; row < nrow; row++) {
1873  double time = gravi_pfits_get_time_acqcam (header, row);
1874  for (int tel = 0; tel < ntel; tel ++)
1875  cpl_table_set (acqcam_table, "TIME", row*ntel+tel, time);
1876  }
1877 
1878  /* Compute mean image */
1879  cpl_image * mean_img = cpl_imagelist_collapse_create (acqcam_imglist);
1880 
1881  /* Compute FIELD columns */
1882  gravi_acqcam_field (mean_img, acqcam_imglist, header,
1883  acqcam_table, o_header, static_param_data);
1884 
1885  CPLCHECK_MSG ("Cannot reduce field images");
1886 
1887 
1888  /* Compute PUPIL columns */
1889  gravi_acqcam_pupil (mean_img, acqcam_imglist, header,
1890  acqcam_table, o_header, static_param_data);
1891  CPLCHECK_MSG ("Cannot reduce pupil images");
1892 
1893  /* Add this output table in the gravi_data */
1894  gravi_data_add_img (output_data, NULL, GRAVI_IMAGING_DATA_ACQ_EXT, mean_img);
1895  CPLCHECK_MSG ("Cannot add acqcam_table in data");
1896 
1897  gravi_data_add_table (output_data, NULL, GRAVI_OI_VIS_ACQ_EXT, acqcam_table);
1898  CPLCHECK_MSG ("Cannot add acqcam_table in data");
1899 
1900  gravi_msg_function_exit(1);
1901  return CPL_ERROR_NONE;
1902 }
1903 
1904 /*----------------------------------------------------------------------------*/
1909 /*----------------------------------------------------------------------------*/
1910 
1911 double gravi_acqcam_z2meter (double PositionPixels, gravi_data * static_param_data)
1912 {
1913  double f_PT ; /* = 14e-3; */ /* pupil tracker lenslet FL*/
1914  double f_lens ; /* = 467e-3; */ /* folding optics lens FL */
1915  double Llambda ; /* = 1.2e-6; */ /* laser diode wavelength */
1916  double D_beam ; /* = 18e-3; */ /* meter */
1917  double D_pixel ; /* = 18e-6; */
1918  double D_AT ; /* = 1.8; */ /* m */
1919  double D_lenslet ; /* = 2 * 1.015e-3; */
1920 
1921  f_PT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS fPT");
1922  cpl_msg_info (cpl_func,"ESO LENS fPT is : %e", f_PT);
1923 
1924  f_lens = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS flens");
1925  cpl_msg_info (cpl_func,"ESO LENS f_lens is : %e", f_lens);
1926 
1927  Llambda = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Llambda");
1928  cpl_msg_info (cpl_func,"ESO LENS Llambda is : %e", Llambda);
1929 
1930  D_beam = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dbeam");
1931  cpl_msg_info (cpl_func,"ESO LENS D_beam is : %e", D_beam);
1932 
1933  D_pixel = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dpixel");
1934  cpl_msg_info (cpl_func,"ESO LENS D_pixel is : %e", D_pixel);
1935 
1936  D_AT = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS DAT");
1937  cpl_msg_info (cpl_func,"ESO LENS D_AT is : %e", D_AT);
1938 
1939  D_lenslet = cpl_propertylist_get_double (gravi_data_get_plist(static_param_data,GRAVI_PRIMARY_HDR_EXT), "ESO LENS Dlenslet");
1940  cpl_msg_info (cpl_func,"ESO LENS D_lenslet is : %e", D_lenslet);
1941 
1942 
1943  double longDef;
1944  longDef = 8 * (f_PT / D_lenslet) * (f_PT / D_lenslet) * 3.5 * D_pixel *
1945  D_beam / (f_PT * D_lenslet) * Llambda / CPL_MATH_2PI * PositionPixels;
1946 
1947  return f_lens * f_lens * longDef / (f_PT + longDef) / f_PT * (D_AT / D_lenslet);
1948 }
1949 
cpl_error_code gravi_acqcam_fit_spot(cpl_image *img, cpl_size ntry, cpl_vector *a, int fitAll, int *nspot)
Fit a pupil spot pattern into an image.
Definition: gravi_acqcam.c:631
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:1443
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, filter median bias.
Definition: gravi_acqcam.c:154
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:1716
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(cpl_propertylist *header, cpl_size tel, cpl_vector *output)
Read the diode position from header into the vector output.
Definition: gravi_acqcam.c:488
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:624
double gravi_image_get_noise_window(cpl_image *img, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
Estimate the median noise in a window. This noise is estimated as the median{img**2}**0.5, so that it is a robust estimate.
Definition: gravi_cpl.c:1608
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
Definition: gravi_data.c:1684
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:1751
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:1972
cpl_error_code gravi_acqcam_get_pup_ref(cpl_propertylist *header, cpl_size tel, cpl_vector *pupref)
Re-arrange the sub-aperture position into the output vector.
Definition: gravi_acqcam.c:550
cpl_error_code gravi_acqcam_spot_imprint(cpl_image *img, cpl_vector *a)
Imprint a cross (pixel=0) in the image of the pupil spot.
Definition: gravi_acqcam.c:433
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:1909
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:925
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:2013
cpl_error_code gravi_reduce_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *static_param_data)
Reduce the ACQ camera images.
double gravi_acqcam_z2meter(double PositionPixels, gravi_data *static_param_data)
Convert z_shift from [pixel] to [meters] Formula extracted from gvacqPupilTracker.c.
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)
Reduce the images of pupil from ACQ.