GRAVI Pipeline Reference Manual 1.8.0
Loading...
Searching...
No Matches
gravity_astrometry.c
Go to the documentation of this file.
1/*
2 * This file is part of the GRAVI Pipeline
3 * Copyright (C) 2002,2003 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 */
19
20#ifdef HAVE_CONFIG_H
21#include <config.h>
22#endif
23
24/*-----------------------------------------------------------------------------
25 Includes
26 -----------------------------------------------------------------------------*/
27
28#include <cpl.h>
29#include <stdio.h>
30#include <string.h>
31#include <time.h>
32
33#include "gravi_data.h"
34#include "gravi_pfits.h"
35#include "gravi_dfs.h"
36#include "gravi_vis.h"
37#include "gravi_utils.h"
38
39#include "gravi_astrometry.h"
40
41/*-----------------------------------------------------------------------------
42 Private function prototypes
43 -----------------------------------------------------------------------------*/
44
45static int gravity_astrometry_create(cpl_plugin *);
46static int gravity_astrometry_exec(cpl_plugin *);
47static int gravity_astrometry_destroy(cpl_plugin *);
48static int gravity_astrometry(cpl_frameset *, cpl_parameterlist *);
49
50/*-----------------------------------------------------------------------------
51 Static variables
52 -----------------------------------------------------------------------------*/
53static char gravity_astro_short[] = "Compute astrometric phase reference";
55 "This recipe computes phase and amplitude referencing for astrometric observations.\n"
56 "It supports on- and off-axis observing strategies, as well as the use of swaps.\n"
58 "* If swaps are present: obtain astrometric solution and compute swap phase reference\n"
59 "* Compute phase reference for the target.\n"
60 "* Write output product with correctly referenced phase.\n"
62 GRAVI_ASTRO_CAL_PHASEREF" : star frames to use for phase referencing\n"
63 GRAVI_ASTRO_TARGET":\tplanet frames to be referenced\n"
64 GRAVI_ASTRO_SWAP":\talternating star/planet frames for swap observing mode\n"
66 GRAVI_ASTRO_PHASE_CALIBRATED" : output astroreduced file with correctly referenced phase\n"
67 "";
68
69/*-----------------------------------------------------------------------------
70 Function code
71 -----------------------------------------------------------------------------*/
72
73/*----------------------------------------------------------------------------*/
83/*----------------------------------------------------------------------------*/
84int cpl_plugin_get_info(cpl_pluginlist * list)
85{
86 cpl_recipe * recipe = cpl_calloc(1, sizeof *recipe );
87 cpl_plugin * plugin = &recipe->interface;
88
89 if (cpl_plugin_init(plugin,
90 CPL_PLUGIN_API,
91 GRAVI_BINARY_VERSION,
92 CPL_PLUGIN_TYPE_RECIPE,
93 "gravity_astrometry",
96 "Calvin Sykes, Mathias Nowak, Sebastian Hoenig",
97 PACKAGE_BUGREPORT,
102 cpl_msg_error(cpl_func, "Plugin initialization failed");
103 (void)cpl_error_set_where(cpl_func);
104 return 1;
105 }
106
107 if (cpl_pluginlist_append(list, plugin)) {
108 cpl_msg_error(cpl_func, "Error adding plugin to list");
109 (void)cpl_error_set_where(cpl_func);
110 return 1;
111 }
112
113 return 0;
114}
115
116/*----------------------------------------------------------------------------*/
124/*----------------------------------------------------------------------------*/
125static int gravity_astrometry_create(cpl_plugin * plugin)
126{
127 cpl_recipe * recipe;
128
129 /* Do not create the recipe if an error code is already set */
130 if (cpl_error_get_code() != CPL_ERROR_NONE) {
131 cpl_msg_error(cpl_func, "%s():%d: An error is already set: %s",
132 cpl_func, __LINE__, cpl_error_get_where());
133 return (int)cpl_error_get_code();
134 }
135
136 if (plugin == NULL) {
137 cpl_msg_error(cpl_func, "Null plugin");
138 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
139 }
140
141 /* Verify plugin type */
142 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
143 cpl_msg_error(cpl_func, "Plugin is not a recipe");
144 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
145 }
146
147 /* Get the recipe */
148 recipe = (cpl_recipe *)plugin;
149
150 /* Create the parameters list in the cpl_recipe object */
151 recipe->parameters = cpl_parameterlist_new();
152 if (recipe->parameters == NULL) {
153 cpl_msg_error(cpl_func, "Parameter list allocation failed");
154 cpl_ensure_code(0, (int)CPL_ERROR_ILLEGAL_OUTPUT);
155 }
156
157 /* Use static names (output_procatg.fits) */
158 gravi_parameter_add_static_name (recipe->parameters);
159
160 /* Astrometry parameters */
161 gravi_parameter_add_astrometry (recipe->parameters);
162
163 return 0;
164}
165
166/*----------------------------------------------------------------------------*/
172/*----------------------------------------------------------------------------*/
173static int gravity_astrometry_exec(cpl_plugin * plugin)
174{
175
176 cpl_recipe * recipe;
177 int recipe_status;
178 cpl_errorstate initial_errorstate = cpl_errorstate_get();
179
180
181 /* Return immediately if an error code is already set */
182 if (cpl_error_get_code() != CPL_ERROR_NONE) {
183 cpl_msg_error(cpl_func, "%s():%d: An error is already set: %s",
184 cpl_func, __LINE__, cpl_error_get_where());
185 return (int)cpl_error_get_code();
186 }
187
188 if (plugin == NULL) {
189 cpl_msg_error(cpl_func, "Null plugin");
190 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
191 }
192
193 /* Verify plugin type */
194 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
195 cpl_msg_error(cpl_func, "Plugin is not a recipe");
196 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
197 }
198
199 /* Get the recipe */
200 recipe = (cpl_recipe *)plugin;
201
202 /* Verify parameter and frame lists */
203 if (recipe->parameters == NULL) {
204 cpl_msg_error(cpl_func, "Recipe invoked with NULL parameter list");
205 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
206 }
207 if (recipe->frames == NULL) {
208 cpl_msg_error(cpl_func, "Recipe invoked with NULL frame set");
209 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
210 }
211
212 /* Invoke the recipe */
213 recipe_status = gravity_astrometry(recipe->frames, recipe->parameters);
214
215 /* Ensure DFS-compliance of the products */
216 if (cpl_dfs_update_product_header(recipe->frames)) {
217 if (!recipe_status){
218 recipe_status = (int)cpl_error_get_code();
219 }
220 }
221
222 if (!cpl_errorstate_is_equal(initial_errorstate)) {
223 /* Dump the error history since recipe execution start.
224 At this point the recipe cannot recover from the error */
225 cpl_errorstate_dump(initial_errorstate, CPL_FALSE, NULL);
226 }
227
228 return recipe_status;
229}
230
231/*----------------------------------------------------------------------------*/
237/*----------------------------------------------------------------------------*/
238static int gravity_astrometry_destroy(cpl_plugin * plugin)
239{
240 cpl_recipe * recipe;
241
242 if (plugin == NULL) {
243 cpl_msg_error(cpl_func, "Null plugin");
244 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
245 }
246
247 /* Verify plugin type */
248 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
249 cpl_msg_error(cpl_func, "Plugin is not a recipe");
250 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
251 }
252
253 /* Get the recipe */
254 recipe = (cpl_recipe *)plugin;
255
256 cpl_parameterlist_delete(recipe->parameters);
257
258 return 0;
259}
260
270static astro_data** load_data(cpl_frameset *frameset, cpl_frameset *used_frameset, gravi_data **out_data)
271{
272 cpl_frame *frame = NULL;
273 gravi_data *tmp_data = NULL;
274 astro_data **data = NULL;
275
276 cpl_size nframes = cpl_frameset_get_size(frameset);
277
278 data = cpl_calloc(nframes, sizeof(astro_data *));
279 for (int i = 0; i < nframes; i++) {
280 frame = cpl_frameset_get_position(frameset, i);
281 tmp_data = gravi_data_load_frame(frame, used_frameset);
282
283 data[i] = gravi_astrometry_load(tmp_data);
284 if (out_data)
285 out_data[i] = tmp_data;
286 else
287 FREE(gravi_data_delete, tmp_data);
288 CPLCHECK_CLEAN("Could not load data");
289 }
290 return data;
291
292cleanup:
293 FREELOOP(gravi_astrometry_delete, data, nframes);
294 return NULL;
295}
296
297/*----------------------------------------------------------------------------*/
304/*----------------------------------------------------------------------------*/
305static int gravity_astrometry(cpl_frameset * frameset,
306 cpl_parameterlist * parlist)
307{
308 cpl_frameset *target_frameset=NULL, *swap_frameset=NULL, *phaseref_frameset=NULL,
309 *used_frameset=NULL;
310
311 cpl_frame *frame=NULL;
312 cpl_size n_target = 0, n_swap = 0, n_phaseref = 0;
313 gravi_data **tgt_data = NULL, **swap_data = NULL;
314 astro_data **tgt_astro = NULL, **swap_astro=NULL, **phaseref_astro=NULL;
315
316 /* Message */
318 cpl_msg_set_time_on();
319 cpl_msg_set_component_on();
321
322 cpl_boolean debug = cpl_parameter_get_bool(
323 cpl_parameterlist_find(parlist, "gravity.astrometry.wait-for-debugger"));
324
325 if (debug) {
326 fprintf(stderr, "PID is: %d\n", getpid());
327 fprintf(stderr, "Waiting for debugger to attach...\n");
328 volatile int attach = 1;
329 while (attach) {
330 if (attach == 0) break;
331 }
332 }
333
334 // const char *target = cpl_parameter_get_string(
335 // cpl_parameterlist_find_const(parlist, "gravity.astrometry.target-name")
336 // );
337 // const char *swap_target = cpl_parameter_get_string(
338 // cpl_parameterlist_find_const(parlist, "gravity.astrometry.swap-target-name")
339 // );
340
341 double ft_mean_flux_threshold = cpl_parameter_get_double(
342 cpl_parameterlist_find(parlist, "gravity.astrometry.ft-mean-flux"));
343
344 cpl_ensure_code(gravi_dfs_set_groups(frameset) == CPL_ERROR_NONE,
345 cpl_error_get_code()) ;
346
347 /* Dispatch the frameset */
348 target_frameset = gravi_frameset_extract_astro_target(frameset);
349 swap_frameset = gravi_frameset_extract_astro_swap(frameset);
350 phaseref_frameset = gravi_frameset_extract_astro_phaseref(frameset);
351
352 // if (cpl_frameset_is_empty(target_frameset)) {
353 // cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
354 // "No ASTRO_TARGET on the frameset");
355 // goto cleanup;
356 // }
357
358 // if ((swap = cpl_frameset_is_empty(swap_frameset))) {
359 // cpl_msg_debug(cpl_func, "No ASTRO_SWAP on the frameset, assuming onaxis");
360 // }
361
362 // if (cpl_frameset_is_empty(phaseref_frameset)) {
363 // cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
364 // "No ASTRO_CAL_PHASEREF on the frameset");
365 // goto cleanup;
366 // }
367
368 /* Insert calibration frame into the used frameset */
369 used_frameset = cpl_frameset_new();
370
371 /* Target data */
372 n_target = cpl_frameset_get_size(target_frameset);
373 tgt_data = cpl_malloc(n_target * sizeof(gravi_data*));
374 tgt_astro = load_data(target_frameset, used_frameset, tgt_data);
375 cpl_msg_debug(cpl_func, "There are %lld ASTRO_TARGET frames", n_target);
376 CPLCHECK_CLEAN("Could not load target data");
377
378 /* Swap data */
379 n_swap = cpl_frameset_get_size(swap_frameset);
380 swap_data = cpl_malloc(n_swap * sizeof(gravi_data*));
381 swap_astro = load_data(swap_frameset, used_frameset, swap_data);
382 cpl_msg_debug(cpl_func, "There are %lld ASTRO_SWAP frames", n_swap);
383 CPLCHECK_CLEAN("Could not load target data");
384
385 /* Data to be used for phase referencing */
386 n_phaseref = cpl_frameset_get_size(phaseref_frameset);
387 phaseref_astro = load_data(phaseref_frameset, used_frameset, NULL);
388 cpl_msg_debug(cpl_func, "There are %lld ASTRO_CAL_PHASEREF frames", n_phaseref);
389 CPLCHECK_CLEAN("Could not load phaseref data");
390
391 /* If there are SWAP files, reduce them first */
392 /* swap files are modified in-place to store the astrometry */
393 if (n_swap > 0)
394 gravi_astrometry_reduce_swaps(swap_astro, n_swap, parlist);
395 CPLCHECK_CLEAN("Could not reduce swaps");
396
397 /* If there are targets, calculate phase reference */
398 /* TODO: if there are no targets, the recipe runs but does nothing */
399 /* TODO: it might be desirable to be able to (f.e.) reduce swaps in isolation */
400 if (n_target > 0) {
401 double ft_mean_flux_tgt = 0.0;
402 for (int i = 0; i < n_target; i++)
403 ft_mean_flux_tgt += gravi_astrometry_get_mean_ftflux(tgt_astro[i]);
404 ft_mean_flux_tgt /= n_target;
405 cpl_msg_debug(cpl_func, "ftOnPlanetMeanFlux=%f", ft_mean_flux_tgt);
406
407 double ft_mean_flux_ref = 0.0;
408 for (int i = 0; i < n_phaseref; i++)
409 ft_mean_flux_ref += gravi_astrometry_get_mean_ftflux(phaseref_astro[i]);
410 ft_mean_flux_ref /= n_phaseref;
411 cpl_msg_debug(cpl_func, "ftOnStarMeanFlux=%f", ft_mean_flux_ref);
412
413 /* Filter based on FT flux threshold, before normalising by coherent FT flux */
414 for (int i = 0; i < n_target; i++) {
415 gravi_astrometry_filter_ftflux(tgt_astro[i], ft_mean_flux_threshold * ft_mean_flux_tgt);
417 }
418
419 for (int i = 0; i < n_phaseref; i++) {
420 gravi_astrometry_filter_ftflux(phaseref_astro[i], ft_mean_flux_threshold * ft_mean_flux_ref);
421 gravi_astrometry_normalise_to_ft(phaseref_astro[i]);
422 }
423
424 cpl_msg_info(cpl_func, "Creating visibility reference from %lld observations", n_phaseref);
425 for (int i = 0; i < n_target; i++)
426 gravi_astrometry_create_phase_reference(tgt_astro[i], phaseref_astro, n_phaseref, swap_astro, n_swap, parlist);
427 CPLCHECK_CLEAN("Could not calculate phase reference");
428
429 for (int i = 0; i < n_target; i++) {
430 cpl_table *phaseref_table = gravi_astrometry_get_phase_reference(tgt_astro[i]);
431
432 frame = cpl_frameset_get_position(target_frameset, i);
433 gravi_data_add_table(tgt_data[i], NULL, "ASTRO_VISREF", phaseref_table);
434 gravi_data_save_new(tgt_data[i], frameset, NULL, NULL, parlist,
435 used_frameset, frame, "gravity_astrometry",
437 CPLCHECK_CLEAN("Could not save ASTRO_VISREF product");
438 }
439 }
440
441 /* Terminate the function */
442 goto cleanup;
443
444cleanup:
445 /* Deallocation of all variables */
446 cpl_msg_info(cpl_func, "Memory cleanup");
447
448 FREE(cpl_frameset_delete, target_frameset);
449 FREE(cpl_frameset_delete, phaseref_frameset);
450 FREE(cpl_frameset_delete, used_frameset);
451
452 FREELOOP(gravi_astrometry_delete, tgt_astro, n_target);
453 FREELOOP(gravi_astrometry_delete, swap_astro, n_swap);
454 FREELOOP(gravi_astrometry_delete, phaseref_astro, n_phaseref);
455
457 return (int)cpl_error_get_code();
458}
void gravi_astrometry_delete(astro_data *self)
cpl_error_code gravi_astrometry_create_phase_reference(astro_data *self, astro_data **phase_refs, cpl_size nphase, astro_data **swaps, cpl_size nswap, cpl_parameterlist *parlist)
Compute the final astrometric phase reference.
cpl_table * gravi_astrometry_get_phase_reference(astro_data *self)
double gravi_astrometry_get_mean_ftflux(astro_data *self)
cpl_error_code gravi_astrometry_filter_ftflux(astro_data *self, double threshold)
Filter based on FT flux threshold and normalise.
cpl_error_code gravi_astrometry_normalise_to_ft(astro_data *self)
Normalise visibilities to average FT flux.
astro_data * gravi_astrometry_load(gravi_data *data)
Load data for astrometry from a gravi_data.
cpl_error_code gravi_astrometry_reduce_swaps(astro_data **swap_data, cpl_size nswap, cpl_parameterlist *parlist)
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
Definition: gravi_data.h:38
#define GRAVI_RECIPE_OUTPUT
Definition: gravi_dfs.h:39
#define GRAVI_ASTRO_TARGET
Definition: gravi_dfs.h:122
#define GRAVI_RECIPE_FLOW
Definition: gravi_dfs.h:37
#define GRAVI_ASTRO_PHASE_CALIBRATED
Definition: gravi_dfs.h:125
#define GRAVI_RECIPE_INPUT
Definition: gravi_dfs.h:38
#define GRAVI_ASTRO_SWAP
Definition: gravi_dfs.h:123
#define GRAVI_ASTRO_CAL_PHASEREF
Definition: gravi_dfs.h:124
cpl_msg_debug(cpl_func, "Spectra has <50 pixels -> don't flat")
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
#define CPLCHECK_CLEAN(msg)
Definition: gravi_utils.h:54
#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 FREELOOP(function, variable, n)
Definition: gravi_utils.h:72
static astro_data ** load_data(cpl_frameset *frameset, cpl_frameset *used_frameset, gravi_data **out_data)
Load input astrometric quantities from ASTROREDUCED file(s).
static char gravity_astro_short[]
int cpl_plugin_get_info(cpl_pluginlist *list)
Build the list of available plugins, for this module.
static int gravity_astrometry_create(cpl_plugin *)
Setup the recipe options
static int gravity_astrometry_destroy(cpl_plugin *)
Destroy what has been created by the 'create' function.
static int gravity_astrometry_exec(cpl_plugin *)
Execute the plugin instance given by the interface.
static int gravity_astrometry(cpl_frameset *, cpl_parameterlist *)
TODO.
static char gravity_astro_description[]
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
gravi_data * gravi_data_load_frame(cpl_frame *frame, cpl_frameset *used_frameset)
Load a FITS file and create a gravi_data.
Definition: gravi_data.c:599
cpl_error_code gravi_data_save_new(gravi_data *self, cpl_frameset *allframes, const char *filename, const char *suffix, const cpl_parameterlist *parlist, cpl_frameset *usedframes, cpl_frame *frame, const char *recipe, cpl_propertylist *applist, const char *proCatg)
Save a gravi data in a CPL-complian FITS file.
Definition: gravi_data.c:925
void gravi_data_delete(gravi_data *self)
Delete a gravi data.
Definition: gravi_data.c:146
cpl_parameter * gravi_parameter_add_static_name(cpl_parameterlist *self)
Definition: gravi_dfs.c:464
cpl_error_code gravi_parameter_add_astrometry(cpl_parameterlist *self)
Definition: gravi_dfs.c:1069
cpl_frameset * gravi_frameset_extract_astro_target(cpl_frameset *frameset)
Definition: gravi_dfs.c:1386
cpl_frameset * gravi_frameset_extract_astro_swap(cpl_frameset *frameset)
Definition: gravi_dfs.c:1390
cpl_error_code gravi_dfs_set_groups(cpl_frameset *set)
Set the group as RAW or CALIB in a frameset.
Definition: gravi_dfs.c:78
cpl_frameset * gravi_frameset_extract_astro_phaseref(cpl_frameset *frameset)
Definition: gravi_dfs.c:1394
void gravity_print_banner(void)
Definition: gravi_dfs.c:61
const char * gravi_get_license(void)
Get the pipeline copyright and license.
Definition: gravi_utils.c:104