65 const char * name_snr,
66 const char * name_gdl);
68 cpl_table * oi_vis2,
const char *name2,
71 const char * name_vis,
72 const char * name_err,
73 const char * name_pha,
74 const char * name_var,
75 const char * name_flag);
77 const char * name_vis,
79 const char * name_flag);
81 const char * name_pha,
82 const char * name_var,
83 const char * name_snr);
85 const char * name_isp,
86 const char * name_gdl,
87 cpl_table * oi_wavelength);
95 double chi2r_threshold,
double chi2r_sigma);
101 cpl_table * vis_ACQ);
104 cpl_table * vis_ACQ);
107 cpl_table * wave_table_sc,
109 cpl_table * wave_table_ft);
114 cpl_table * wave_table_sc,
115 cpl_table * wave_table_ft,
116 cpl_propertylist *
header,
117 const cpl_parameterlist * parlist);
121 cpl_table * wave_table,
122 cpl_table * disp_table,
123 cpl_propertylist *
header,
124 const cpl_parameterlist * parlist);
127 cpl_table * wave_table,
128 cpl_propertylist *
header,
129 const cpl_parameterlist * parlist);
151 const char * name_snr,
152 const char * name_gdl)
155 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
156 cpl_ensure_code (name_snr, CPL_ERROR_NULL_INPUT);
157 cpl_ensure_code (name_gdl, CPL_ERROR_NULL_INPUT);
160 double * snr = cpl_table_get_data_double (oi_vis, name_snr);
161 double * gdl = cpl_table_get_data_double (oi_vis, name_gdl);
166 cpl_size nrow = cpl_table_get_nrow (oi_vis) / 6;
167 for (
int base=0; base<6; base++ ) {
170 for (
int tri=0; tri<2; tri++ ) {
178 cpl_msg_debug(cpl_func,
"Found triangle tels %i%i -> bases (%i,%i, %i,%i)",
182 for (cpl_size row = 0; row < nrow; row++) {
184 double snrN = CPL_MIN( snr[row*6+b1], snr[row*6+b2] );
185 if ( snrN > snr[row*6+base] ) {
186 snr[row*6+base] = snrN;
187 gdl[row*6+base] = sign1 * gdl[row*6+b1] + sign2 * gdl[row*6+b2];
194 cpl_msg_debug(cpl_func,
"TODO: loop on 3-baseline bootstrap");
199 return CPL_ERROR_NONE;
219 cpl_table * oi_vis2,
const char *name2,
223 cpl_ensure_code (oi_vis1, CPL_ERROR_NULL_INPUT);
224 cpl_ensure_code (oi_vis2, CPL_ERROR_NULL_INPUT);
225 cpl_ensure_code (name1, CPL_ERROR_NULL_INPUT);
226 cpl_ensure_code (name2, CPL_ERROR_NULL_INPUT);
227 cpl_ensure_code (phasediff, CPL_ERROR_ILLEGAL_OUTPUT);
231 cpl_type type1 = cpl_table_get_column_type (oi_vis1, name1 );
232 cpl_type type2 = cpl_table_get_column_type (oi_vis2, name2 );
233 cpl_size nrow1 = cpl_table_get_nrow (oi_vis1) / nbase;
234 cpl_size nrow2 = cpl_table_get_nrow (oi_vis2) / nbase;
236 if ( type1 != CPL_TYPE_DOUBLE_COMPLEX ||
237 type2 != CPL_TYPE_DOUBLE_COMPLEX ||
239 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
"input columns not conformables or not DOUBLE COMPLEX");
241 double complex * data1 = cpl_table_get_data_double_complex (oi_vis1, name1);
242 double complex * data2 = cpl_table_get_data_double_complex (oi_vis2, name2);
246 for (
int base = 0; base < nbase; base++) {
247 double complex phasor = 0.0;
250 for (cpl_size row = 0; row < nrow1; row++) {
251 phasor += data1[row*nbase+base] * conj( data2[row*6+base] );
254 phasediff[base] = carg( phasor );
257 for (cpl_size row = 0; row < nrow1; row++) {
258 data1[row*6+base] *= cexp( -I*phasediff[base] );
263 return CPL_ERROR_NONE;
281 const char * name_vis,
282 const char * name_err,
283 const char * name_pha,
284 const char * name_var,
285 const char * name_flag)
288 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
289 cpl_ensure_code (name_vis, CPL_ERROR_NULL_INPUT);
290 cpl_ensure_code (name_err, CPL_ERROR_NULL_INPUT);
291 cpl_ensure_code (name_pha, CPL_ERROR_ILLEGAL_OUTPUT);
292 cpl_ensure_code (name_var, CPL_ERROR_ILLEGAL_OUTPUT);
293 cpl_ensure_code (name_flag,CPL_ERROR_ILLEGAL_OUTPUT);
295 cpl_size nrow = cpl_table_get_nrow (oi_vis);
299 double complex * phasorRaw = cpl_table_get_data_double_complex (oi_vis, name_pha);
302 double * varRaw = cpl_table_get_data_double (oi_vis, name_var);
305 cpl_array** tVis = cpl_table_get_data_array (oi_vis, name_vis);
306 cpl_array** tErr = cpl_table_get_data_array (oi_vis, name_err);
307 cpl_array** tFlag = cpl_table_get_data_array (oi_vis, name_flag);
308 cpl_ensure_code (tFlag, CPL_ERROR_ILLEGAL_INPUT);
311 int nwave = cpl_array_get_size (tVis[0]);
316 for (cpl_size row = 0; row < nrow; row++ ) {
318 double complex * cpx = cpl_array_get_data_double_complex (tVis[row]);
319 double complex * err = cpl_array_get_data_double_complex (tErr[row]);
320 cpl_ensure_code (cpx, CPL_ERROR_ILLEGAL_INPUT);
321 cpl_ensure_code (err, CPL_ERROR_ILLEGAL_INPUT);
322 phasorRaw[row] = 0.0;
325 int * flag = cpl_array_get_data_int (tFlag[row]);
329 for (cpl_size wave = 0; wave < nwave; wave++)
332 phasorRaw[row] += cpx[wave];
333 varRaw[row] += cabs(err[wave]) * cabs(err[wave]);
342 return CPL_ERROR_NONE;
360 const char * name_vis,
361 const char * name_is,
362 const char * name_flag)
365 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
366 cpl_ensure_code (name_vis, CPL_ERROR_NULL_INPUT);
367 cpl_ensure_code (name_is, CPL_ERROR_NULL_INPUT);
368 cpl_ensure_code (name_flag, CPL_ERROR_NULL_INPUT);
370 cpl_size nrow = cpl_table_get_nrow (oi_vis);
373 cpl_array** tVis = cpl_table_get_data_array (oi_vis, name_vis);
374 cpl_ensure_code (tVis, CPL_ERROR_ILLEGAL_INPUT);
376 cpl_array** tFlag = cpl_table_get_data_array (oi_vis, name_flag);
377 cpl_ensure_code (tFlag, CPL_ERROR_ILLEGAL_INPUT);
380 int nwave = cpl_array_get_size (tVis[0]);
386 double complex * interSpectreRaw = cpl_table_get_data_double_complex (oi_vis, name_is);
391 for (cpl_size row = 0; row < nrow; row++) {
393 double complex * cpx = cpl_array_get_data_double_complex (tVis[row]);
394 cpl_ensure_code (cpx, CPL_ERROR_ILLEGAL_INPUT);
396 int * flag = cpl_array_get_data_int (tFlag[row]);
397 cpl_ensure_code (flag, CPL_ERROR_ILLEGAL_INPUT);
399 interSpectreRaw[row] = 0.0;
403 for (cpl_size wave = 0; wave < nwave-1; wave++)
405 if (!flag[wave] && !flag[wave+1]) {
406 interSpectreRaw[row] += cpx[wave] * conj( cpx[wave+1] );
415 return CPL_ERROR_NONE;
433 const char * name_pha,
434 const char * name_var,
435 const char * name_snr)
438 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
439 cpl_ensure_code (name_pha, CPL_ERROR_NULL_INPUT);
440 cpl_ensure_code (name_var, CPL_ERROR_NULL_INPUT);
441 cpl_ensure_code (name_snr, CPL_ERROR_ILLEGAL_OUTPUT);
443 cpl_size nrow = cpl_table_get_nrow (oi_vis);
446 double complex * cpx = cpl_table_get_data_double_complex (oi_vis, name_pha);
447 double * dbl = cpl_table_get_data_double (oi_vis, name_var);
449 CPLCHECK_MSG (
"Cannot load data to compute SNR column");
453 double * snr = cpl_table_get_data_double (oi_vis, name_snr);
458 for (cpl_size row = 0; row < nrow; row++) {
459 snr[row] = cabs( cpx[row] ) / sqrt (fabs(dbl[row]));
465 return CPL_ERROR_NONE;
483 const char * name_isp,
484 const char * name_gdl,
485 cpl_table * oi_wavelength)
488 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
489 cpl_ensure_code (name_isp, CPL_ERROR_NULL_INPUT);
490 cpl_ensure_code (name_gdl, CPL_ERROR_NULL_INPUT);
491 cpl_ensure_code (oi_wavelength, CPL_ERROR_NULL_INPUT);
493 cpl_size nrow = cpl_table_get_nrow (oi_vis);
496 double complex * cpx1 = cpl_table_get_data_double_complex (oi_vis, name_isp);
498 CPLCHECK_MSG(
"Cannot load data to compute GDELAY column");
502 double * gdl = cpl_table_get_data_double (oi_vis, name_gdl);
506 cpl_size wave = cpl_table_get_ncol (oi_wavelength)/2;
507 double factor = 1./fabs(1./cpl_table_get (oi_wavelength,
"EFF_WAVE", wave, NULL) -
508 1./cpl_table_get (oi_wavelength,
"EFF_WAVE", wave+1, NULL));
510 cpl_msg_debug (cpl_func,
"Compute the coherence length (%e m)", factor);
515 for (cpl_size row = 0; row < nrow; row++ ) {
517 gdl[row] = (double)carg( cpx1[row] ) / (2.0*M_PI) * factor;
523 return CPL_ERROR_NONE;
537 const cpl_parameterlist * parlist)
540 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
541 cpl_ensure_code (parlist, CPL_ERROR_NULL_INPUT);
547 cpl_ensure_code (p2vmred_header, CPL_ERROR_ILLEGAL_INPUT);
550 for (
int type_data = 0; type_data < 2; type_data ++) {
553 for (
int pol = 0; pol < npol; pol++) {
557 cpl_size nrow = cpl_table_get_nrow (oi_flux) /
ntel;
558 cpl_size nwave = cpl_table_get_column_depth (oi_flux,
"FLUX");
576 for (cpl_size wave = 0; wave < nwave ; wave ++) {
577 double value = (double)cpl_array_get (array, wave, NULL) / nrow;
578 if (value > 0) stat[0]++;
579 if (value > 0.5) stat[1]++;
582 FREE (cpl_array_delete, array);
586 sprintf (qc_name,
"ESO QC OUTLIER_SC_P%i", pol+1);
587 cpl_propertylist_update_int (p2vmred_header, qc_name, stat[0]);
588 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"Channels with at least one outlier");
591 sprintf (qc_name,
"ESO QC OUTLIER50_SC_P%i",pol+1);
592 cpl_propertylist_update_int (p2vmred_header, qc_name, stat[1]);
593 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"Channels with more than 50% outlier");
603 return CPL_ERROR_NONE;
621 const cpl_parameterlist * parlist)
625 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
641 for (
int type_data = 0; type_data < 2; type_data ++) {
644 cpl_msg_info (cpl_func,
"Cannot compute real-time snr for %s "
658 nsmooth = 1./ periode_sc;
659 nsmooth = CPL_MIN (CPL_MAX (nsmooth, 0), 20);
664 for (
int pol = 0; pol < npol; pol++) {
666 cpl_msg_info (cpl_func,
"Compute SNR for pol %i of %s (smoothed over %i frames)",
697 for (
int type_data = 0; type_data < 2; type_data ++) {
700 cpl_msg_info (cpl_func,
"Cannot compute bootstraped snr %s "
707 cpl_msg_info(cpl_func,
"Compute bootstraped GDELAY_BOOT and SNR for %s...",(type_data==
GRAVI_FT?
"FT":
"SC"));
711 cpl_table * oi_vis_p2, * oi_vis_p1;
717 cpl_table_duplicate_column (oi_vis_p1,
"IPHASOR_BOOT", oi_vis_p1,
"IPHASOR_SMT");
721 cpl_table_duplicate_column (oi_vis_p1,
"PHASOR_BOOT", oi_vis_p1,
"PHASOR_SMT");
724 cpl_table_duplicate_column (oi_vis_p1,
"PHASOR_VAR_BOOT", oi_vis_p1,
"PHASOR_VAR_SMT");
730 cpl_msg_info(cpl_func,
"Add the signal of both polarisation to enhance SNR and GDELAY accuracy");
750 cpl_msg_info(cpl_func,
"Add the phase difference as QC parameters");
751 for (
int base = 0; base < nbase; base++) {
752 sprintf (qc_name,
"ESO QC PHASE_POLDIFF_%s%d%d",(type_data==
GRAVI_FT?
"FT":
"SC"),
755 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[deg] differential phase" );
756 cpl_msg_info (cpl_func,
"%s=%f [deg]",qc_name,phasediff[base] * CPL_MATH_DEG_RAD);
760 CPLCHECK_MSG(
"Cannot add columns of second polarisation");
770 int nstabilize = (type_data ==
GRAVI_FT ? 50 : 1);
774 CPLCHECK_MSG(
"Cannot compute and fill SNR_BOOT or GDELAY_BOOT column");
778 cpl_msg_info (cpl_func,
"Duplicate columns in polarisation 2");
779 if (cpl_table_has_column (oi_vis_p2,
"SNR_BOOT"))
780 cpl_table_erase_column(oi_vis_p2,
"SNR_BOOT");
781 cpl_table_duplicate_column (oi_vis_p2,
"SNR_BOOT", oi_vis_p1,
"SNR_BOOT");
783 if (cpl_table_has_column (oi_vis_p2,
"GDELAY_BOOT"))
784 cpl_table_erase_column(oi_vis_p2,
"GDELAY_BOOT");
785 cpl_table_duplicate_column (oi_vis_p2,
"GDELAY_BOOT", oi_vis_p1,
"GDELAY_BOOT");
786 CPLCHECK_MSG(
"Cannot duplicate column for polarisation 2");
791 cpl_msg_info (cpl_func,
"Erase useless columns in P2VMREDUCED");
793 cpl_table_erase_column (oi_vis_p1,
"IPHASOR");
794 cpl_table_erase_column (oi_vis_p1,
"IPHASOR_SMT");
795 cpl_table_erase_column (oi_vis_p1,
"IPHASOR_BOOT");
796 cpl_table_erase_column (oi_vis_p1,
"PHASOR");
797 cpl_table_erase_column (oi_vis_p1,
"PHASOR_SMT");
798 cpl_table_erase_column (oi_vis_p1,
"PHASOR_BOOT");
799 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR");
800 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR_SMT");
801 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR_BOOT");
802 cpl_table_erase_column (oi_vis_p1,
"SNR_BOOT_TMP");
806 cpl_table_erase_column (oi_vis_p2,
"IPHASOR");
807 cpl_table_erase_column (oi_vis_p2,
"IPHASOR_SMT");
808 cpl_table_erase_column (oi_vis_p2,
"PHASOR");
809 cpl_table_erase_column (oi_vis_p2,
"PHASOR_SMT");
810 cpl_table_erase_column (oi_vis_p2,
"PHASOR_VAR");
811 cpl_table_erase_column (oi_vis_p2,
"PHASOR_VAR_SMT");
816 for (cpl_size base = 0; base < nbase; base++) {
820 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean bootstrapped SNR");
826 return CPL_ERROR_NONE;
850 cpl_table * vis_FT,
int nbase_ft,
854 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
855 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
858 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase_sc;
859 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase_ft;
865 sprintf (full_name,
"FIRST_%s", name);
867 int * first_ft = cpl_table_get_data_int (vis_SC, full_name);
869 sprintf (full_name,
"LAST_%s", name);
871 int * last_ft = cpl_table_get_data_int (vis_SC, full_name);
873 sprintf (full_name,
"NFRAME_%s", name);
875 int * nframe_ft = cpl_table_get_data_int (vis_SC, full_name);
880 for (cpl_size base_sc = 0; base_sc < nbase_sc; base_sc++) {
885 cpl_size base_ft = base_sc % nbase_ft;
888 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
891 while ( cpl_table_get (vis_FT,
"TIME", row_ft * nbase_ft + base_ft, NULL) <
892 cpl_table_get (vis_SC,
"TIME", row_sc * nbase_sc + base_sc, NULL) - dit_sc/2.) {
893 if (row_ft >= nrow_ft-1)
break;
896 first_ft[row_sc * nbase_sc + base_sc] = row_ft;
899 while ( cpl_table_get (vis_FT,
"TIME", row_ft * nbase_ft + base_ft, NULL) <
900 cpl_table_get (vis_SC,
"TIME", row_sc * nbase_sc + base_sc, NULL) + dit_sc/2.) {
901 if (row_ft >= nrow_ft-1)
break;
904 last_ft[row_sc * nbase_sc + base_sc] = row_ft;
907 if ( first_ft[row_sc * nbase_sc + base_sc] < 2 ||
908 last_ft[row_sc * nbase_sc + base_sc] > nrow_ft - 2 ) {
909 if (!strcmp (name,
"ACQ"))
911 cpl_msg_info (cpl_func,
"Not enough %s data to synchronise with DIT %lli over %lli", name, row_sc+1, nrow_sc);
913 cpl_msg_warning (cpl_func,
"Not enough %s data to synchronise with DIT %lli over %lli", name, row_sc+1, nrow_sc);
915 first_ft[row_sc * nbase_sc + base_sc] = 0;
916 last_ft[row_sc * nbase_sc + base_sc] = 0;
920 nframe_ft[row_sc * nbase_sc + base_sc] = last_ft[row_sc * nbase_sc + base_sc] - first_ft[row_sc * nbase_sc + base_sc];
926 return CPL_ERROR_NONE;
946 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
947 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
950 cpl_size nbase = 6,
ntel = 4, nclo = 4;
951 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
955 double * pFactor = cpl_table_get_data_double (vis_SC,
"P_FACTOR");
958 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
959 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
962 double * flux = cpl_table_get_data_double (flux_FT,
"TOTALFLUX_UNSMOOTHED");
967 for (cpl_size base = 0; base < nbase; base++) {
968 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
969 cpl_size nsc = row_sc * nbase + base;
972 double sf0f1 = 0.0, f0 = 0.0, f1 = 0.0;
973 for (cpl_size row_ft = first_ft[nsc]; row_ft < last_ft[nsc]; row_ft++) {
976 if (flux[t0] < 0.0 || flux[t1] < 0.0)
978 sf0f1 += sqrt (flux[t0] * flux[t1]);
984 if (f0==0 || f1==0)
continue;
987 pFactor[nsc] = sf0f1 * sf0f1 / (f0 * f1);
993 double * p3Factor = cpl_table_get_data_double (t3_SC,
"P3_FACTOR");
996 first_ft = cpl_table_get_data_int (t3_SC,
"FIRST_FT");
997 last_ft = cpl_table_get_data_int (t3_SC,
"LAST_FT");
1000 for (cpl_size closure = 0; closure < nclo; closure++) {
1001 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
1002 cpl_size nt3 = row_sc * nclo + closure;
1005 double sf0f1f2 = 0.0, f0 = 0.0, f1 = 0.0, f2 = 0.0;
1006 for (cpl_size row_ft = first_ft[nt3]; row_ft < last_ft[nt3]; row_ft++) {
1010 if (flux[t0] < 0.0 || flux[t1] < 0.0 || flux[t2] < 0.0)
1012 sf0f1f2 += pow (flux[t0] * flux[t1] * flux[t2], 1.0 / 3);
1019 if (f0==0 || f1==0 || f2==0)
continue;
1022 p3Factor[nt3] = sf0f1f2 * sf0f1f2 * sf0f1f2 / (f0 * f1 * f2);
1027 return CPL_ERROR_NONE;
1046 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1047 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
1050 cpl_size nbase = 6,
ntel = 4;
1051 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1054 cpl_size nwave_ft = cpl_table_get_column_depth (flux_FT,
"FLUX");
1055 cpl_array **p_flux = cpl_table_get_data_array (flux_FT,
"FLUX");
1059 cpl_array **p_pFactor = cpl_table_get_data_array (vis_FT,
"P_FACTOR");
1063 double *sf0f1 = cpl_malloc(nwave_ft *
sizeof(
double));
1064 double *f0 = cpl_malloc(nwave_ft *
sizeof(
double));
1065 double *f1 = cpl_malloc(nwave_ft *
sizeof(
double));
1068 for (cpl_size base = 0; base < nbase; base++) {
1069 for (cpl_size row_ft = 0; row_ft < nrow_ft; row_ft++) {
1070 cpl_size nft = row_ft * nbase + base;
1073 int window_start = CPL_MAX(0, row_ft - window_width);
1074 int window_end = CPL_MIN(row_ft + window_width + 1, nrow_ft);
1076 double *pFactor = cpl_array_get_data_double(p_pFactor[nft]);
1078 for (cpl_size wave = 0; wave < nwave_ft; wave++)
1079 sf0f1[wave] = f0[wave] = f1[wave] = 0.0;
1081 for (cpl_size row_ft = window_start; row_ft < window_end; row_ft++) {
1085 double *pf0 = cpl_array_get_data_double(p_flux[t0]);
1086 double *pf1 = cpl_array_get_data_double(p_flux[t1]);
1089 for (cpl_size wave = 0; wave < nwave_ft; wave++) {
1090 if (pf0[wave] < 0.0 || pf1[wave] < 0.0)
1093 sf0f1[wave] += sqrt (pf0[wave] * pf1[wave]);
1094 f0[wave] += pf0[wave];
1095 f1[wave] += pf1[wave];
1099 for (cpl_size wave = 0; wave < nwave_ft; wave++) {
1100 if (f0[wave] == 0 || f1[wave] == 0)
1104 pFactor[wave] = sf0f1[wave] * sf0f1[wave] / (f0[wave] * f1[wave]);
1114 return CPL_ERROR_NONE;
1133 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1134 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1136 cpl_size nbase = 6,
ntel = 4;
1137 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) /
ntel;
1138 cpl_size nwave_sc = cpl_table_get_column_depth (flux_SC,
"FLUX");
1141 cpl_array ** flux_sc = cpl_table_get_data_array (flux_SC,
"FLUX");
1147 cpl_array ** f1f2_sc = cpl_table_get_data_array (vis_SC,
"F1F2");
1149 for (cpl_size base = 0; base < nbase; base++) {
1150 for (cpl_size row = 0; row < nrow_sc; row ++) {
1154 f1f2_sc[row*nbase+base] = cpl_array_cast (flux_sc[row*
ntel+
GRAVI_BASE_TEL[base][0]], CPL_TYPE_DOUBLE);
1155 cpl_array_multiply (f1f2_sc[row*nbase+base], flux_sc[row*
ntel+
GRAVI_BASE_TEL[base][1]]);
1160 return CPL_ERROR_NONE;
1181 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1182 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
1184 cpl_size
ntel = 4, nbase = 6;
1185 cpl_size nrow_ft = cpl_table_get_nrow (flux_FT) /
ntel;
1186 cpl_size nwave_ft = cpl_table_get_column_depth (flux_FT,
"FLUX");
1189 double * total_flux_ft = cpl_table_get_data_double (flux_FT,
"TOTALFLUX");
1190 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
1195 cpl_array ** mean_spectra = cpl_malloc (4 *
sizeof (cpl_array*));
1197 for (cpl_size tel = 0; tel <
ntel; tel ++) {
1198 mean_spectra[tel] = cpl_array_duplicate (flux_ft[tel]);
1199 for (cpl_size n = 1; n < nrow_ft; n ++) cpl_array_add (mean_spectra[tel], flux_ft[n*
ntel+tel]);
1200 cpl_array_divide_scalar (mean_spectra[tel], cpl_array_get_mean (mean_spectra[tel]) * nwave_ft);
1207 cpl_array ** f1f2_ft = cpl_table_get_data_array (vis_FT,
"F1F2");
1212 for (cpl_size base = 0; base < nbase; base++) {
1215 for (cpl_size n = 0; n < nrow_ft; n ++) {
1216 f1f2_ft[n*nbase+base] = cpl_array_duplicate (mean_spectra[t0]);
1217 cpl_array_multiply (f1f2_ft[n*nbase+base], mean_spectra[t1]);
1218 cpl_array_multiply_scalar (f1f2_ft[n*nbase+base], total_flux_ft[n*
ntel+t0] * total_flux_ft[n*
ntel+t1]);
1224 FREELOOP (cpl_array_delete, mean_spectra, 4);
1227 return CPL_ERROR_NONE;
1245 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1249 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1250 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
1251 cpl_array ** visData_ft = cpl_table_get_data_array (vis_FT,
"VISDATA");
1257 cpl_array ** phaseref_ft = cpl_table_get_data_array (vis_FT,
"SELF_REF");
1263 for (cpl_size base = 0; base < nbase; base ++) {
1266 for (cpl_size n = 3; n < nrow_ft - 3; n ++) {
1267 phaseref_ft[n*nbase+base] = cpl_array_cast (visData_ft[(n-1)*nbase+base], CPL_TYPE_DOUBLE_COMPLEX);
1268 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+1)*nbase+base]);
1269 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], 2.0);
1270 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n-2)*nbase+base]);
1271 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+2)*nbase+base]);
1272 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], 2.0);
1273 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n-3)*nbase+base]);
1274 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+3)*nbase+base]);
1276 cpl_array_arg (phaseref_ft[n*nbase+base]);
1277 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], -1.0);
1279 CPLCHECK_MSG(
"Cannot compute the PHASE_REF for the FT");
1281 phaseref_ft[0*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1282 phaseref_ft[1*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1283 phaseref_ft[2*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1284 phaseref_ft[(nrow_ft-3)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1285 phaseref_ft[(nrow_ft-2)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1286 phaseref_ft[(nrow_ft-1)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1290 cpl_msg_warning (cpl_func,
"Change of PHASE_REF_FT common channel new... to be decided !!");
1293 return CPL_ERROR_NONE;
1310 cpl_table * vis_ACQ)
1313 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1314 cpl_ensure_code (vis_ACQ, CPL_ERROR_NULL_INPUT);
1316 cpl_size nbase = 6,
ntel = 4;
1317 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1320 int * first = cpl_table_get_data_int (vis_SC,
"FIRST_ACQ");
1321 int * last = cpl_table_get_data_int (vis_SC,
"LAST_ACQ");
1325 int * pup_n = cpl_table_get_data_int (vis_ACQ,
"PUPIL_NSPOT");
1326 double * pup_x = cpl_table_get_data_double (vis_ACQ,
"PUPIL_X");
1327 double * pup_y = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Y");
1328 double * pup_z = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Z");
1333 double * pup_x_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_X");
1335 double * pup_y_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_Y");
1337 double * pup_z_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_Z");
1340 for (cpl_size base = 0; base < nbase; base++) {
1341 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1342 cpl_size nsc = row_sc * nbase + base;
1345 cpl_size nframe = 0;
1346 for (cpl_size row = first[nsc] ; row < last[nsc]; row++) {
1350 if (pup_n[row0] != 0 && pup_n[row1] !=0 ) {
1351 pup_x_sc[nsc] += pup_x[row0] - pup_x[row1];
1352 pup_y_sc[nsc] += pup_y[row0] - pup_y[row1];
1353 pup_z_sc[nsc] += pup_z[row0] - pup_z[row1];
1362 pup_x_sc[nsc] /= (double)nframe;
1363 pup_y_sc[nsc] /= (double)nframe;
1364 pup_z_sc[nsc] /= (double)nframe;
1371 return CPL_ERROR_NONE;
1388 cpl_table * vis_ACQ)
1391 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1392 cpl_ensure_code (vis_ACQ, CPL_ERROR_NULL_INPUT);
1395 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) /
ntel;
1398 int * first = cpl_table_get_data_int (flux_SC,
"FIRST_ACQ");
1399 int * last = cpl_table_get_data_int (flux_SC,
"LAST_ACQ");
1403 int * pup_n = cpl_table_get_data_int (vis_ACQ,
"PUPIL_NSPOT");
1404 double * pup_x = cpl_table_get_data_double (vis_ACQ,
"PUPIL_X");
1405 double * pup_y = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Y");
1406 double * pup_z = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Z");
1411 double * pup_x_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_X");
1413 double * pup_y_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_Y");
1415 double * pup_z_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_Z");
1418 for (cpl_size tel = 0; tel <
ntel; tel++) {
1419 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1420 cpl_size nsc = row_sc *
ntel + tel;
1423 cpl_size nframe = 0;
1424 for (cpl_size row = first[nsc] ; row < last[nsc]; row++) {
1425 cpl_size row0 = row *
ntel + tel;
1427 if (pup_n[row0] != 0) {
1428 pup_x_sc[nsc] += pup_x[row0];
1429 pup_y_sc[nsc] += pup_y[row0];
1430 pup_z_sc[nsc] += pup_z[row0];
1439 pup_x_sc[nsc] /= (double)nframe;
1440 pup_y_sc[nsc] /= (double)nframe;
1441 pup_z_sc[nsc] /= (double)nframe;
1448 return CPL_ERROR_NONE;
1465 cpl_table * wave_table)
1468 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1469 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1470 cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
1473 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1477 cpl_size nwave_sc = cpl_table_get_column_depth (vis_SC,
"VISDATA");
1478 double * twopi_wavenumber_sc = cpl_malloc (nwave_sc *
sizeof(
double));
1479 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
1480 twopi_wavenumber_sc[wave] = CPL_MATH_2PI / cpl_table_get (wave_table,
"EFF_WAVE", wave, NULL);
1484 int * first_met = cpl_table_get_data_int (vis_SC,
"FIRST_MET");
1485 int * last_met = cpl_table_get_data_int (vis_SC,
"LAST_MET");
1490 double * opd_met_fc = cpl_table_get_data_double (vis_MET,
"OPD_FC");
1491 cpl_array ** opd_met_tel = cpl_table_get_data_array (vis_MET,
"OPD_TEL");
1493 double * opd_met_fc_corr = cpl_table_get_data_double (vis_MET,
"OPD_FC_CORR");
1494 double * opd_met_telfc_mcorr = cpl_table_get_data_double (vis_MET,
"OPD_TELFC_MCORR");
1495 cpl_array ** opd_met_telfc_corr = cpl_table_get_data_array (vis_MET,
"OPD_TELFC_CORR");
1501 cpl_array ** phase_metdit_telfc = cpl_table_get_data_array (vis_SC,
"PHASE_MET_TELFC");
1504 double * opd_metdit_fc = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC");
1507 cpl_array ** opd_metdit_tel = cpl_table_get_data_array (vis_SC,
"OPD_MET_TEL");
1510 double * opd_metdit_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
1513 double * opd_metdit_telfc_mcorr = cpl_table_get_data_double (vis_SC,
"OPD_MET_TELFC_MCORR");
1516 cpl_array ** opd_metdit_telfc_corr = cpl_table_get_data_array (vis_SC,
"OPD_MET_TELFC_CORR");
1521 cpl_msg_info (cpl_func,
"Compute OPD_MET_TEL, OPD_MET_FC_CORR, OPD_MET_TELFC_MCORR...");
1523 for (cpl_size base = 0; base < nbase; base++) {
1524 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1525 cpl_size nsc = row_sc * nbase + base;
1531 _Complex
double * phasor_metdit_telfc_ptr = cpl_array_get_data_double_complex(phasor_metdit_telfc);
1534 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1539 double opd_astro = opd_met_fc_corr[nmet0] - opd_met_fc_corr[nmet1] + opd_met_telfc_mcorr[nmet0] - opd_met_telfc_mcorr[nmet1];
1541 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
1542 phasor_metdit_telfc_ptr[wave] = phasor_metdit_telfc_ptr[wave] + cexp(opd_astro * I * twopi_wavenumber_sc[wave]);
1546 opd_metdit_fc_corr[nsc] += opd_met_fc_corr[nmet0] - opd_met_fc_corr[nmet1];
1547 opd_metdit_telfc_mcorr[nsc] += opd_met_telfc_mcorr[nmet0] - opd_met_telfc_mcorr[nmet1];
1550 cpl_array_add (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet0]);
1551 cpl_array_subtract (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet1]);
1554 cpl_array_add (opd_metdit_tel[nsc], opd_met_tel[nmet0]);
1555 cpl_array_subtract (opd_metdit_tel[nsc], opd_met_tel[nmet1]);
1558 opd_metdit_fc[nsc] += opd_met_fc[nmet0] - opd_met_fc[nmet1];
1564 cpl_size nframe = last_met[nsc] - first_met[nsc];
1566 opd_metdit_fc_corr[nsc] /= nframe;
1567 opd_metdit_telfc_mcorr[nsc] /= nframe;
1568 cpl_array_divide_scalar (opd_metdit_telfc_corr[nsc], (
double)nframe);
1569 cpl_array_divide_scalar (opd_metdit_tel[nsc], (
double)nframe);
1570 opd_metdit_fc[nsc] /= nframe;
1573 cpl_array_arg (phasor_metdit_telfc);
1574 phase_metdit_telfc[nsc] = cpl_array_cast(phasor_metdit_telfc, CPL_TYPE_DOUBLE);
1576 cpl_array_delete(phasor_metdit_telfc);
1577 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1583 if (cpl_table_has_column (vis_MET,
"FIELD_FIBER_DX")) {
1585 cpl_msg_info (cpl_func,
"Compute the information coming from VIS_ACQ");
1587 double * fdx_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DX");
1588 double * fdy_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DY");
1591 double * fdx_metdit = cpl_table_get_data_double (vis_SC,
"FIELD_FIBER_DX");
1594 double * fdy_metdit = cpl_table_get_data_double (vis_SC,
"FIELD_FIBER_DY");
1597 for (cpl_size base = 0; base < nbase; base++) {
1598 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1599 cpl_size nsc = row_sc * nbase + base;
1602 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1607 fdx_metdit[nsc] += fdx_met[nmet0] - fdx_met[nmet1];
1608 fdy_metdit[nsc] += fdy_met[nmet0] - fdy_met[nmet1];
1610 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1613 cpl_size nframe = last_met[nsc] - first_met[nsc];
1615 fdx_metdit[nsc] /= nframe;
1616 fdy_metdit[nsc] /= nframe;
1624 if (cpl_table_has_column (vis_MET,
"OPD_PUPIL") && cpl_table_has_column (vis_MET,
"OPD_TTPUP")) {
1626 cpl_msg_info (cpl_func,
"Compute mean OPD_MET_PUPIL, OPD_MET_PUPIL_STDDEV and OPD_MET_TTPUP...");
1628 double * opd_met_pupil = cpl_table_get_data_double (vis_MET,
"OPD_PUPIL");
1629 double * opd_met_ttpup = cpl_table_get_data_double (vis_MET,
"OPD_TTPUP");
1632 double * opd_metdit_pupil = cpl_table_get_data_double (vis_SC,
"OPD_MET_PUPIL");
1635 double * opd_metdit_pupil_stddev = cpl_table_get_data_double (vis_SC,
"OPD_MET_PUPIL_STDDEV");
1638 double * opd_metdit_ttpup = cpl_table_get_data_double (vis_SC,
"OPD_MET_TTPUP");
1643 for (cpl_size base = 0; base < nbase; base++) {
1644 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1645 cpl_size nsc = row_sc * nbase + base;
1648 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1653 opd_metdit_pupil[nsc] += opd_met_pupil[nmet0] - opd_met_pupil[nmet1];
1656 opd_metdit_ttpup[nsc] += opd_met_ttpup[nmet0] - opd_met_ttpup[nmet1];
1662 cpl_size nframe = last_met[nsc] - first_met[nsc];
1664 opd_metdit_pupil[nsc] /= nframe;
1665 opd_metdit_ttpup[nsc] /= nframe;
1669 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1674 double tmp = opd_met_pupil[nmet0] - opd_met_pupil[nmet1] - opd_metdit_pupil[nsc];
1675 opd_metdit_pupil_stddev[nsc] += tmp * tmp;
1682 opd_metdit_pupil_stddev[nsc] = sqrt(opd_metdit_pupil_stddev[nsc] / (nframe-1));
1685 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1690 FREE (cpl_free, twopi_wavenumber_sc);
1692 return CPL_ERROR_NONE;
1720 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1721 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1723 cpl_size nbase = 6,
ntel = 4;
1724 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1725 cpl_size nrow_met = cpl_table_get_nrow (vis_MET) /
ntel;
1730 double periode_ft = cpl_table_get (vis_FT,
"TIME", nbase, NULL) -
1731 cpl_table_get (vis_FT,
"TIME", 0, NULL);
1732 cpl_msg_info (cpl_func,
"PERIOD FT = %g [us]", periode_ft);
1735 double * phase_met_fc = cpl_table_get_data_double (vis_MET,
"PHASE_FC_DRS");
1741 double * phase_metdit_fc = cpl_table_get_data_double (vis_FT,
"PHASE_MET_FC");
1746 for (cpl_size base = 0; base < nbase; base++) {
1750 for (cpl_size last_met = 0, first_met = 0, row_ft = 0; row_ft < nrow_ft; row_ft ++) {
1751 cpl_size nft = row_ft * nbase + base;
1752 double time_ft = cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL);
1761 first_met = CPL_MAX (CPL_MIN (last_met - 5, nrow_met - 1), 0);
1762 while ((cpl_table_get (vis_MET,
"TIME", first_met*
ntel, NULL) <= (time_ft + 0.0))) {
1764 if (first_met == nrow_met)
break;
1769 last_met = CPL_MAX (CPL_MIN (first_met - 1, nrow_met - 1), 0);
1770 while ((cpl_table_get (vis_MET,
"TIME", last_met*
ntel, NULL) < (time_ft + periode_ft))) {
1772 if (last_met == nrow_met)
break;
1777 if (row_ft < 5 && last_met == 0) last_met = 1;
1781 if (row_ft > nrow_ft-5 && first_met == nrow_met) first_met = nrow_met - 1;
1782 if (row_ft > nrow_ft-5 && last_met == nrow_met) last_met = nrow_met;
1785 if ( last_met == 0 || last_met > nrow_met ) {
1786 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
1787 "Not enough MET data to synchronise "
1788 "with FT DIT %lli over %lli", row_ft+1, nrow_ft);
1791 if ( last_met - first_met > 0 ) {
1795 for (cpl_size rin_met = first_met; rin_met < last_met; rin_met ++)
1796 phase_metdit_fc[nft] += phase_met_fc[rin_met*
ntel+tel0] - phase_met_fc[rin_met*
ntel+tel1];
1798 phase_metdit_fc[nft] /= (last_met - first_met);
1802 cpl_size rowa_met = first_met-1;
1803 cpl_size rowb_met = last_met;
1805 double phia_met = phase_met_fc[rowa_met*
ntel+tel0] - phase_met_fc[rowa_met*
ntel+tel1];
1806 double phib_met = phase_met_fc[rowb_met*
ntel+tel0] - phase_met_fc[rowb_met*
ntel+tel1];
1807 double timea_met = cpl_table_get (vis_MET,
"TIME",rowa_met*
ntel, NULL);
1808 double timeb_met = cpl_table_get (vis_MET,
"TIME",rowb_met*
ntel, NULL);
1810 phase_metdit_fc[nft] = phia_met +
1811 (phib_met - phia_met) * (time_ft - timea_met) / (timeb_met - timea_met);
1815 CPLCHECK_MSG (
"Fail to compute metrology per FT base from metrology per tel");
1821 return CPL_ERROR_NONE;
1848 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1849 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1852 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1853 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1858 double * phase_sc = cpl_table_get_data_double (vis_SC,
"OPD");
1864 double * phase_scdit = cpl_table_get_data_double (vis_FT,
"OPD_SC");
1869 for (cpl_size base = 0; base < nbase; base++) {
1870 cpl_size row_ft = 0;
1872 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1873 cpl_size nsc = row_sc * nbase + base;
1874 double time_sc = cpl_table_get (vis_SC,
"TIME", nsc, NULL);
1877 while ((cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL) < (time_sc - dit_sc/2))) {
1879 if (row_ft >= nrow_ft)
break;
1883 while ((cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL) < (time_sc + dit_sc/2))) {
1885 phase_scdit[row_ft*nbase+base] = phase_sc[nsc];
1888 if (row_ft >= nrow_ft)
break;
1891 if (row_ft >= nrow_ft) {
1892 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
1893 "Not enough FT data to synchronise "
1894 "with SC DIT %lli over %lli", row_sc+1, nrow_sc);
1902 return CPL_ERROR_NONE;
1922 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1923 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1926 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) /
ntel;
1929 int * first_met = cpl_table_get_data_int (flux_SC,
"FIRST_MET");
1930 int * last_met = cpl_table_get_data_int (flux_SC,
"LAST_MET");
1935 double * opd_met_fc = cpl_table_get_data_double (vis_MET,
"OPD_FC");
1936 cpl_array ** opd_met_tel = cpl_table_get_data_array (vis_MET,
"OPD_TEL");
1938 double * opd_met_fc_corr = cpl_table_get_data_double (vis_MET,
"OPD_FC_CORR");
1939 double * opd_met_telfc_mcorr = cpl_table_get_data_double (vis_MET,
"OPD_TELFC_MCORR");
1940 cpl_array ** opd_met_telfc_corr = cpl_table_get_data_array (vis_MET,
"OPD_TELFC_CORR");
1941 double * phase_met_fc = cpl_table_get_data_double (vis_MET,
"PHASE_FC_DRS");
1942 double * phase_met_fcft = cpl_table_get_data_double (vis_MET,
"PHASE_FCFT_DRS");
1943 double * phase_met_fcsc = cpl_table_get_data_double (vis_MET,
"PHASE_FCSC_DRS");
1950 double * phase_metdit_fc = cpl_table_get_data_double (flux_SC,
"PHASE_MET_FC");
1953 double * phase_metdit_fcft = cpl_table_get_data_double (flux_SC,
"PHASE_MET_FCFT");
1956 double * phase_metdit_fcsc = cpl_table_get_data_double (flux_SC,
"PHASE_MET_FCSC");
1959 double * opd_metdit_fc = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC");
1962 cpl_array ** opd_metdit_tel = cpl_table_get_data_array (flux_SC,
"OPD_MET_TEL");
1965 double * opd_metdit_fc_corr = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC_CORR");
1968 double * opd_metdit_telfc_mcorr = cpl_table_get_data_double (flux_SC,
"OPD_MET_TELFC_MCORR");
1971 cpl_array ** opd_metdit_telfc_corr = cpl_table_get_data_array (flux_SC,
"OPD_MET_TELFC_CORR");
1976 for (cpl_size tel = 0; tel <
ntel; tel++) {
1977 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1978 cpl_size nsc = row_sc *
ntel + tel;
1984 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1985 cpl_size nmet = row_met *
ntel + tel;
1988 cpl_array_add (opd_metdit_tel[nsc], opd_met_tel[nmet]);
1991 opd_metdit_fc[nsc] += opd_met_fc[nmet];
1994 opd_metdit_fc_corr[nsc] += opd_met_fc_corr[nmet];
1995 opd_metdit_telfc_mcorr[nsc] += opd_met_telfc_mcorr[nmet];
1998 cpl_array_add (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet]);
2001 phase_metdit_fc[nsc] += phase_met_fc[nmet];
2004 phase_metdit_fcft[nsc] += phase_met_fcft[nmet];
2007 phase_metdit_fcsc[nsc] += phase_met_fcsc[nmet];
2013 cpl_size nframe = last_met[nsc] - first_met[nsc];
2015 opd_metdit_fc_corr[nsc] /= nframe;
2016 opd_metdit_telfc_mcorr[nsc] /= nframe;
2017 phase_metdit_fc[nsc] /= nframe;
2018 phase_metdit_fcft[nsc] /= nframe;
2019 phase_metdit_fcsc[nsc] /= nframe;
2020 cpl_array_divide_scalar (opd_metdit_telfc_corr[nsc], (
double)nframe);
2022 cpl_array_divide_scalar (opd_metdit_tel[nsc], (
double)nframe);
2023 opd_metdit_fc[nsc] /= nframe;
2033 if (cpl_table_has_column (vis_MET,
"FIELD_FIBER_DX")) {
2035 double * fdx_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DX");
2036 double * fdy_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DY");
2039 double * fdx_metdit = cpl_table_get_data_double (flux_SC,
"FIELD_FIBER_DX");
2042 double * fdy_metdit = cpl_table_get_data_double (flux_SC,
"FIELD_FIBER_DY");
2045 for (cpl_size tel = 0; tel <
ntel; tel++) {
2046 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2047 cpl_size nsc = row_sc *
ntel + tel;
2050 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
2051 cpl_size nmet0 = row_met *
ntel + tel;
2054 fdx_metdit[nsc] += fdx_met[nmet0];
2055 fdy_metdit[nsc] += fdy_met[nmet0];
2057 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
2060 cpl_size nframe = last_met[nsc] - first_met[nsc];
2062 fdx_metdit[nsc] /= nframe;
2063 fdy_metdit[nsc] /= nframe;
2071 return CPL_ERROR_NONE;
2091 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
2092 cpl_ensure_code (fddl_table, CPL_ERROR_NULL_INPUT);
2095 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) /
ntel;
2098 int * first_fddl = cpl_table_get_data_int (flux_SC,
"FIRST_FDDL");
2099 int * last_fddl = cpl_table_get_data_int (flux_SC,
"LAST_FDDL");
2101 cpl_array ** ftpos = cpl_table_get_data_array (fddl_table,
"FT_POS");
2102 cpl_array ** scpos = cpl_table_get_data_array (fddl_table,
"SC_POS");
2103 cpl_array ** oplair = cpl_table_get_data_array (fddl_table,
"OPL_AIR");
2109 double * ftpos_flux = cpl_table_get_data_double (flux_SC,
"FT_POS");
2112 double * scpos_flux = cpl_table_get_data_double (flux_SC,
"SC_POS");
2115 double * oplair_flux = cpl_table_get_data_double (flux_SC,
"OPL_AIR");
2120 for (cpl_size tel = 0; tel <
ntel; tel++) {
2121 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2122 cpl_size nsc = row_sc *
ntel + tel;
2125 for (cpl_size row_fddl = first_fddl[nsc] ; row_fddl < last_fddl[nsc]; row_fddl++) {
2126 ftpos_flux[nsc] += cpl_array_get (ftpos[row_fddl], tel, NULL);
2127 scpos_flux[nsc] += cpl_array_get (scpos[row_fddl], tel, NULL);
2128 oplair_flux[nsc] += cpl_array_get (oplair[row_fddl], tel, NULL);
2132 cpl_size nframe = last_fddl[nsc] - first_fddl[nsc];
2134 ftpos_flux[nsc] /= nframe;
2135 scpos_flux[nsc] /= nframe;
2136 oplair_flux[nsc] /= nframe;
2143 return CPL_ERROR_NONE;
2162 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
2163 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
2166 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) /
ntel;
2167 cpl_size nwave_sc = cpl_table_get_column_depth (flux_SC,
"FLUX");
2168 cpl_size nwave_ft = cpl_table_get_column_depth (flux_FT,
"FLUX");
2171 int * first_ft = cpl_table_get_data_int (flux_SC,
"FIRST_FT");
2172 int * last_ft = cpl_table_get_data_int (flux_SC,
"LAST_FT");
2174 cpl_array ** flag_sc = cpl_table_get_data_array (flux_SC,
"FLAG");
2175 cpl_array ** flux_sc = cpl_table_get_data_array (flux_SC,
"FLUX");
2176 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
2182 double * total_flux_scdit = cpl_table_get_data_double (flux_SC,
"TOTALFLUX_SC");
2185 double * total_flux_ftdit = cpl_table_get_data_double (flux_SC,
"TOTALFLUX_FT");
2189 for (cpl_size tel = 0; tel <
ntel; tel++) {
2190 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2191 cpl_size nsc = row_sc *
ntel + tel;
2194 total_flux_scdit[nsc] = 0.0;
2196 for (
int wave = 0; wave < nwave_sc; wave++) {
2197 if (!cpl_array_get (flag_sc[nsc], wave, NULL)) {
2198 total_flux_scdit[nsc] += cpl_array_get (flux_sc[nsc], wave, NULL);
2205 total_flux_scdit[nsc] *= (double)nwave_sc / (
double)nvalid;
2209 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
2210 total_flux_ftdit[nsc] += cpl_array_get_mean (flux_ft[row_ft *
ntel + tel]) * nwave_ft;
2213 CPLCHECK_MSG(
"Issue in the loop to average FT and SC flux per frame");
2219 return CPL_ERROR_NONE;
2241 cpl_table * wave_table_sc,
2243 cpl_table * wave_table_ft)
2246 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2247 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
2248 cpl_ensure_code (wave_table_sc, CPL_ERROR_NULL_INPUT);
2249 cpl_ensure_code (wave_table_ft, CPL_ERROR_NULL_INPUT);
2253 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2254 cpl_size nwave_sc = cpl_table_get_column_depth (vis_SC,
"VISDATA");
2255 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
2258 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
2259 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
2261 cpl_array ** visErr_ft = cpl_table_get_data_array (vis_FT,
"VISERR");
2262 cpl_array ** visData_ft = cpl_table_get_data_array (vis_FT,
"VISDATA");
2267 double meanwave_ft = cpl_table_get_column_mean (wave_table_ft,
"EFF_WAVE");
2268 double * wavenumber_sc = cpl_malloc (nwave_sc *
sizeof(
double));
2269 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
2270 wavenumber_sc[wave] = 1. / cpl_table_get (wave_table_sc,
"EFF_WAVE", wave, NULL);
2277 cpl_array ** visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA_FT");
2280 cpl_array ** visVar_ftdit = cpl_table_get_data_array (vis_SC,
"VISVAR_FT");
2283 cpl_array ** visPower_ftdit = cpl_table_get_data_array (vis_SC,
"VISPOWER_FT");
2286 cpl_array ** vFactor = cpl_table_get_data_array (vis_SC,
"V_FACTOR");
2289 cpl_array ** vFactor_ftdit = cpl_table_get_data_array (vis_SC,
"V_FACTOR_FT");
2292 double * vFactor_wl = cpl_table_get_data_double (vis_SC,
"V_FACTOR_WL");
2297 for (cpl_size base = 0; base < nbase; base++) {
2298 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
2299 int nsc = row_sc * nbase + base;
2306 double vFactor_incoh = 0.0;
2307 double complex vFactor_coh = 0.0 + I * 0.0;
2308 double vFactor_var = 0.0;
2311 cpl_size nframe = last_ft[nsc] - first_ft[nsc];
2312 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
2313 cpl_array * tmp_data;
2316 tmp_data = visErr_ft[row_ft * nbase + base];
2317 double *visVar_ftdit_ptr = cpl_array_get_data_double(visVar_ftdit[nsc]);
2318 _Complex
double *tmp_data_ptr = cpl_array_get_data_double_complex(tmp_data);
2319 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
2320 visVar_ftdit_ptr[wave] = visVar_ftdit_ptr[wave] +
gravi_pow2(cabs(tmp_data_ptr[wave]));
2324 tmp_data = visData_ft[row_ft * nbase + base];
2325 _Complex
double *visData_ftdit_ptr = cpl_array_get_data_double_complex(visData_ftdit[nsc]);
2326 tmp_data_ptr = cpl_array_get_data_double_complex(tmp_data);
2328 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
2329 visData_ftdit_ptr[wave] = visData_ftdit_ptr[wave] + tmp_data_ptr[wave];
2333 double *visPower_ftdit_ptr = cpl_array_get_data_double(visPower_ftdit[nsc]);
2334 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
2335 visPower_ftdit_ptr[wave] = visPower_ftdit_ptr[wave] +
gravi_pow2(cabs(tmp_data_ptr[wave]));
2340 vFactor_incoh +=
gravi_pow2(cabs (cpl_array_get_mean_complex (visData_ft[row_ft * nbase + base])) * nwave_ft);
2341 vFactor_coh += cpl_array_get_mean_complex (visData_ft[row_ft * 6 + base]) * nwave_ft;
2342 _Complex
double *visErr_ft_ptr = cpl_array_get_data_double_complex(visErr_ft[row_ft * nbase + base]);
2343 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2344 vFactor_var +=
gravi_pow2(cabs (visErr_ft_ptr[ wave]));
2347 CPLCHECK_MSG(
"Issue in the loop to build average FT quantities");
2352 vFactor_ftdit[nsc] = cpl_array_new (nwave_ft, CPL_TYPE_DOUBLE);
2353 _Complex
double *visData_ftdit_ptr = cpl_array_get_data_double_complex(visData_ftdit[nsc]);
2355 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2356 cpl_array_set (vFactor_ftdit[nsc], wave,
2357 (
gravi_pow2(cabs (visData_ftdit_ptr[ wave])) -
2358 cpl_array_get (visVar_ftdit[nsc], wave, NULL)) /
2359 (cpl_array_get (visPower_ftdit[nsc], wave, NULL) -
2360 cpl_array_get (visVar_ftdit[nsc], wave, NULL)) / (
double)nframe);
2365 vFactor_wl[nsc] = (cabs(vFactor_coh)*cabs(vFactor_coh) - vFactor_var) /
2366 (vFactor_incoh - vFactor_var) / (double)nframe;
2370 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2371 cpl_array_set (vFactor_ftdit[nsc], wave,0);
2373 vFactor_wl[nsc] = 0;
2378 vFactor[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2379 double a2 = log (CPL_MAX(CPL_MIN (vFactor_wl[nsc], 1.0),1e-10)) * meanwave_ft*meanwave_ft;
2380 for (cpl_size wave = 0; wave < nwave_sc; wave ++) {
2381 cpl_array_set (vFactor[nsc], wave, exp ( a2 *
gravi_pow2(wavenumber_sc[wave])));
2387 FREE (cpl_free, wavenumber_sc);
2390 return CPL_ERROR_NONE;
2413 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2414 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
2418 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2421 int * reject_flag_ft = cpl_table_get_data_int (vis_FT,
"REJECTION_FLAG");
2422 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
2423 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
2429 double * fringedet_ftdit = cpl_table_get_data_double (vis_SC,
"FRINGEDET_RATIO");
2434 for (cpl_size base = 0; base < nbase; base++) {
2435 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
2436 int nsc = row_sc * nbase + base;
2439 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
2440 fringedet_ftdit[nsc] += (reject_flag_ft[row_ft * nbase + base] == 0 ? 1 : 0);
2444 cpl_size nframe = last_ft[nsc] - first_ft[nsc];
2445 if (nframe != 0) fringedet_ftdit[nsc] /= (double)nframe;
2450 return CPL_ERROR_NONE;
2473 cpl_table * wavesc_table,
2474 cpl_table * waveft_table,
2475 cpl_propertylist *
header,
2476 const cpl_parameterlist * parlist)
2479 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2480 cpl_ensure_code (wavesc_table, CPL_ERROR_NULL_INPUT);
2482 cpl_array ** visData_ftdit;
2483 double * gdelay_ftdit;
2484 const char * output_name = NULL, * coeff_name = NULL;
2485 if (waveft_table != NULL) {
2486 cpl_msg_info (cpl_func,
"Compute reference phase of SC from the FT data");
2487 visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA_FT");
2488 gdelay_ftdit = cpl_table_get_data_double (vis_SC,
"GDELAY_FT");
2489 output_name =
"PHASE_REF";
2490 coeff_name =
"PHASE_REF_COEFF";
2492 cpl_msg_info (cpl_func,
"Compute reference phase of SC from the SC");
2493 visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA");
2494 gdelay_ftdit = cpl_table_get_data_double (vis_SC,
"GDELAY");
2495 waveft_table = wavesc_table;
2496 output_name =
"SELF_REF";
2497 coeff_name =
"SELF_REF_COEFF";
2502 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2503 cpl_size nwave_sc = cpl_table_get_nrow (wavesc_table);
2504 cpl_size nwave_ft = cpl_table_get_nrow (waveft_table);
2511 cpl_size mindeg = 0;
2519 cpl_msg_info (cpl_func,
"phaseref with polynomial mindeg=%lli to maxdeg=%lli", mindeg, maxdeg);
2520 cpl_propertylist_update_int (
header,
"ESO QC PHASEREF_SC MINDEG", mindeg);
2521 cpl_propertylist_set_comment (
header,
"ESO QC PHASEREF_SC MINDEG",
"fit of FT phase");
2522 cpl_propertylist_update_int (
header,
"ESO QC PHASEREF_SC MAXDEG", maxdeg);
2523 cpl_propertylist_set_comment (
header,
"ESO QC PHASEREF_SC MAXDEG",
"fit of FT phase");
2526 cpl_polynomial * fit = cpl_polynomial_new (1);
2529 cpl_matrix * sigma_ft = cpl_matrix_new (1,nwave_ft);
2530 cpl_vector * wave_sc = cpl_vector_new (nwave_sc);
2531 cpl_array * wavenumber_ft = cpl_array_new (nwave_ft, CPL_TYPE_DOUBLE);
2533 double lbd0 = cpl_table_get_column_mean (wavesc_table,
"EFF_WAVE");
2534 double delta0 = cpl_table_get_column_max (wavesc_table,
"EFF_WAVE") -
2535 cpl_table_get_column_min (wavesc_table,
"EFF_WAVE");
2536 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2537 double lbd = cpl_table_get (waveft_table,
"EFF_WAVE", wave, NULL);
2538 cpl_matrix_set (sigma_ft, 0, wave, (lbd0/lbd - 1.) * lbd0/delta0 );
2539 cpl_array_set (wavenumber_ft, wave, 1./lbd);
2541 for (cpl_size wave = 0; wave < nwave_sc; wave ++) {
2542 cpl_vector_set (wave_sc, wave, cpl_table_get (wavesc_table,
"EFF_WAVE", wave, NULL));
2549 cpl_array ** phaseref = cpl_table_get_data_array (vis_SC, output_name);
2553 cpl_array ** phase_coeff = cpl_table_get_data_array (vis_SC, coeff_name);
2558 for (cpl_size base = 0; base < nbase; base++) {
2559 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2560 int nsc = row_sc * nbase + base;
2564 cpl_array * phaseref_ftdit = cpl_array_cast (visData_ftdit[nsc], CPL_TYPE_DOUBLE_COMPLEX);
2568 double mean_phase = carg (cpl_array_get_mean_complex (phaseref_ftdit));
2569 cpl_array_multiply_scalar_complex (phaseref_ftdit, cexp(- I * mean_phase));
2572 cpl_array_arg (phaseref_ftdit);
2574 cpl_array_add_scalar (phaseref_ftdit, mean_phase);
2580 cpl_vector * input = cpl_vector_wrap (nwave_ft, cpl_array_get_data_double (phaseref_ftdit));
2581 cpl_polynomial_fit (fit, sigma_ft, NULL, input, NULL, CPL_FALSE, &mindeg, &maxdeg);
2582 cpl_vector_unwrap (input);
2583 cpl_array_delete (phaseref_ftdit);
2586 phase_coeff[nsc] = cpl_array_new (maxdeg+1, CPL_TYPE_DOUBLE);
2587 for (cpl_size d = 0; d < maxdeg+1; d++)
2588 cpl_array_set (phase_coeff[nsc], d, cpl_polynomial_get_coeff (fit, &d));
2591 phaseref[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2592 for (cpl_size w = 0; w < nwave_sc; w++) {
2593 double delta = (lbd0/cpl_vector_get(wave_sc, w) - 1.) * lbd0/delta0 ;
2594 cpl_array_set (phaseref[nsc], w, cpl_polynomial_eval_1d (fit, delta, NULL));
2598 cpl_array_multiply_scalar (phaseref[nsc], -1.0);
2604 FREE (cpl_vector_delete, wave_sc);
2605 FREE (cpl_matrix_delete, sigma_ft);
2606 FREE (cpl_array_delete, wavenumber_ft);
2607 FREE (cpl_polynomial_delete, fit);
2610 return CPL_ERROR_NONE;
2630 double chi2r_threshold,
2634 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
2635 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2637 cpl_msg_info (cpl_func,
"chi2_threshold = %f", chi2r_threshold);
2638 cpl_msg_info (cpl_func,
"chi2_sigma = %f", chi2r_sigma);
2641 cpl_size nrow = cpl_table_get_nrow (flux_SC) /
ntel;
2642 cpl_size nwave = cpl_table_get_column_depth (flux_SC,
"FLUX");
2644 cpl_array ** chi2 = cpl_table_get_data_array (flux_SC,
"CHI2");
2647 cpl_array ** flux_outliers = cpl_table_get_data_array (flux_SC,
"FLAG");
2652 cpl_vector * vector = cpl_vector_new (nwave);
2653 cpl_vector * med = NULL;
2656 for (cpl_size row = 0; row < nrow; row ++) {
2659 for (cpl_size wave = 0; wave < nwave ; wave ++) {
2660 double value = cpl_array_get (chi2[row*
ntel+0], wave, NULL);
2661 cpl_vector_set (vector, wave, value);
2667 cpl_vector_add_scalar(med, cpl_vector_get_mean(med)*1e-9);
2668 cpl_vector_divide (vector, med);
2669 FREE (cpl_vector_delete, med);
2672 for (cpl_size wave = 0; wave < nwave ; wave ++) {
2673 double value = cpl_vector_get (vector, wave);
2674 if (value > chi2r_threshold) cpl_array_set (flux_outliers[row*
ntel+0], wave, 1);
2679 cpl_vector_subtract_scalar (vector, 1.0);
2683 cpl_vector_add_scalar(med, cpl_vector_get_mean(med)*1e-9);
2684 cpl_vector_divide (vector, med);
2685 FREE (cpl_vector_delete, med);
2688 for (cpl_size wave = 0; wave < nwave ; wave ++) {
2689 double value = cpl_vector_get (vector, wave);
2690 if (value > chi2r_sigma) cpl_array_set (flux_outliers[row*
ntel+0], wave, 1);
2695 FREE (cpl_vector_delete, vector);
2703 for (cpl_size row = 0; row < nrow; row ++) {
2704 for (cpl_size tel = 1; tel <
ntel; tel++) {
2705 cpl_array_add (flux_outliers[row*
ntel+tel], flux_outliers[row*
ntel+0]);
2715 cpl_array ** vis_outliers = cpl_table_get_data_array (vis_SC,
"FLAG");
2719 for (cpl_size row = 0; row < nrow; row ++) {
2720 for (cpl_size base = 0; base < nbase; base++) {
2721 cpl_array_add (vis_outliers[row*nbase+base], flux_outliers[row*
ntel+0]);
2728 return CPL_ERROR_NONE;
2745 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
2746 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
2751 return CPL_ERROR_NONE;
2766 cpl_table * disp_table)
2769 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
2772 cpl_size nrow = cpl_table_get_nrow (flux_SC) /
ntel;
2783 if (disp_table == NULL) {
2785 return CPL_ERROR_NONE;
2789 cpl_msg_info (cpl_func,
"Load the linearity model coeficients");
2793 cpl_size disp_order = cpl_table_get_column_depth (disp_table,
"LIN_FDDL_SC");
2797 double * ftpos = cpl_table_get_data_double (flux_SC,
"FT_POS");
2798 double * scpos = cpl_table_get_data_double (flux_SC,
"SC_POS");
2799 double * fddl = cpl_table_get_data_double (flux_SC,
"FDDL");
2800 double * fddl_ft = cpl_table_get_data_double (flux_SC,
"FDDL_FT");
2801 double * fddl_sc = cpl_table_get_data_double (flux_SC,
"FDDL_SC");
2805 for (cpl_size tel = 0; tel <
ntel; tel++) {
2806 for (cpl_size row = 0; row < nrow; row ++) {
2807 cpl_size nsc = row *
ntel + tel;
2813 for (
int o = 0; o < disp_order; o++) {
2814 fddl_sc[nsc] += lin_fddl_sc[tel][o] * pow (scpos[nsc], (
double)o) * 1.0e-6;
2815 fddl_ft[nsc] += lin_fddl_ft[tel][o] * pow (ftpos[nsc], (
double)o) * 1.0e-6;
2817 fddl[nsc]=0.5*(fddl_ft[nsc]+fddl_sc[nsc]);
2822 FREE (cpl_free, lin_fddl_sc);
2823 FREE (cpl_free, lin_fddl_ft);
2826 return CPL_ERROR_NONE;
2846 cpl_table * flux_SC,
2847 cpl_table * wave_table,
2848 cpl_table * disp_table,
2849 cpl_propertylist *
header,
2850 const cpl_parameterlist * parlist)
2853 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2854 cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
2855 cpl_ensure_code (
header, CPL_ERROR_NULL_INPUT);
2857 cpl_size nwave_sc = cpl_table_get_column_depth (vis_SC,
"VISDATA");
2860 if (disp_table == NULL) {
2861 cpl_msg_info (cpl_func,
"Cannot create OPD_DISP, not DISP_MODEL table");
2863 return CPL_ERROR_NONE;
2869 cpl_msg_info (cpl_func,
"Load the dispersion model coeficients");
2872 double ** n_mean = cpl_malloc (4 *
sizeof(
double*));
2873 double ** n_diff = cpl_malloc (4 *
sizeof(
double*));
2874 for (
int t = 0; t < 4; t++) {
2875 n_mean[t] = cpl_calloc (nwave_sc,
sizeof(
double));
2876 n_diff[t] = cpl_calloc (nwave_sc,
sizeof(
double));
2879 if ( cpl_table_has_column (disp_table,
"N_MEAN") &&
2880 cpl_table_has_column (disp_table,
"WAVE0")) {
2883 cpl_msg_info (cpl_func,
"Dispersion model as N/Nmet [wave0/wave-1]");
2884 cpl_size disp_order = cpl_table_get_column_depth (disp_table,
"N_MEAN");
2885 for (
int t = 0; t < 4; t++) {
2886 double lbd0 = cpl_table_get (disp_table,
"WAVE0", t, NULL);
2887 for (
int w = 0; w < nwave_sc; w++ ) {
2888 double lbd = cpl_table_get (wave_table,
"EFF_WAVE", w, NULL);
2889 double xfit = lbd0/lbd - 1.0;
2890 for (
int o = 0; o < disp_order; o++) {
2899 cpl_msg_error (cpl_func,
"The DISP_MODEL is not recognized... contact the DRS team");
2900 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
2901 "The DISP_MODEL is not recognized... contact the DRS team");
2906 cpl_size nbase = 6,
ntel = 4;
2907 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2911 cpl_array ** opd_disp = cpl_table_get_data_array (vis_SC,
"OPD_DISP");
2914 double * gd_disp = cpl_table_get_data_double (vis_SC,
"GDELAY_DISP");
2917 cpl_array ** phase_disp = cpl_table_get_data_array (vis_SC,
"PHASE_DISP");
2922 double * opd_met = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC");
2923 double * fddl = cpl_table_get_data_double (flux_SC,
"FDDL");
2928 double * wavenumber_sc = cpl_malloc (nwave_sc *
sizeof(
double));
2929 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
2930 wavenumber_sc[wave] = 1. / cpl_table_get (wave_table,
"EFF_WAVE", wave, NULL);
2936 double * opl_zero_fc = cpl_malloc (4 *
sizeof(
double));
2940 for (t = 0; t <
ntel; t++) {
2941 opl_zero_fc[t] = 0.0;
2944 cpl_msg_info (cpl_func,
"Metrology zero calculation is enabled!");
2947 for (t = 0; t <
ntel; t++) {
2948 sprintf (name,
"ESO OCS MET OPL_ZERO_FC%i", t+1);
2949 if (cpl_propertylist_has (
header, name)) {
2950 opl_zero_fc[t] = cpl_propertylist_get_double (
header, name)*1e-3;
2951 sprintf (name,
"ESO FDDL MET OFFSET%i", t+1);
2952 opl_zero_fc[t] += cpl_propertylist_get_double (
header, name)*1e-3;
2953 cpl_msg_info (cpl_func,
"Updating metrology zero with OCS MET OPL_ZERO_FC%i and FDDL MET OFFSET%i: %f [mm]", t+1, t+1, opl_zero_fc[t]);
2958 for (t = 0; t <
ntel; t++) {
2962 * (wavenumber_sc[nwave_sc/2-1] - wavenumber_sc[nwave_sc/2+1])
2963 / (n_mean[t][nwave_sc/2-1] * wavenumber_sc[nwave_sc/2-1] - n_mean[t][nwave_sc/2+1] * wavenumber_sc[nwave_sc/2+1]);
2964 cpl_msg_info (cpl_func,
"Updating metrology zero with QC/PRO MET GD_ZERO_FC%i: %f [mm]", t+1, opl_zero_fc[t]);
2968 for (t = 0; t <
ntel; t++) {
2971 cpl_msg_info (cpl_func,
"Updating metrology zero with QC/PRO MET OPL_ZERO_FC%i: %f [mm]", t+1, opl_zero_fc[t]);
2975 cpl_msg_info (cpl_func,
"Metrology zero calculation is disabled!");
2980 for (cpl_size base = 0; base < nbase; base++) {
2981 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2983 cpl_size nsc = row_sc * nbase + base;
2989 cpl_size wave = nwave_sc / 2;
2990 double s1 = wavenumber_sc[wave-1];
2991 double s2 = wavenumber_sc[wave+1];
2993 (n_mean[t0][wave-1] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][wave-1] * fddl[t0f]) -
2994 (n_mean[t1][wave-1] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][wave-1] * fddl[t1f]);
2996 (n_mean[t0][wave+1] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][wave+1] * fddl[t0f]) -
2997 (n_mean[t1][wave+1] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][wave+1] * fddl[t1f]);
2999 gd_disp[nsc] = (o1*s1 - o2*s2) / (s1-s2);
3003 phase_disp[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
3005 for (
int w = 0; w < nwave_sc ; w++) {
3006 cpl_array_set (phase_disp[nsc], w,
3007 ((n_mean[t0][w] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][w] * fddl[t0f]) -
3008 (n_mean[t1][w] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][w] * fddl[t1f]) -
3009 gd_disp[nsc]) * wavenumber_sc[w] * CPL_MATH_2PI);
3017 opd_disp[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
3019 for (
int w = 0; w < nwave_sc ; w++) {
3020 cpl_array_set (opd_disp[nsc], w,
3021 (n_mean[t0][w] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][w] * fddl[t0f]) -
3022 (n_mean[t1][w] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][w] * fddl[t1f]) );
3031 FREE (cpl_free, wavenumber_sc);
3032 FREE (cpl_free, opl_zero_fc);
3035 return CPL_ERROR_NONE;
3053 cpl_table * wave_table,
3054 cpl_propertylist *
header,
3055 const cpl_parameterlist * parlist)
3058 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
3059 cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
3060 cpl_ensure_code (
header, CPL_ERROR_NULL_INPUT);
3063 cpl_size nwave = cpl_table_get_column_depth (vis_SC,
"VISDATA");
3064 cpl_size nrow = cpl_table_get_nrow (vis_SC);
3067 double * ucoord = cpl_table_get_data_double (vis_SC,
"UCOORD");
3068 double * vcoord = cpl_table_get_data_double (vis_SC,
"VCOORD");
3069 cpl_array ** phase_ref = cpl_table_get_data_array (vis_SC,
"PHASE_REF");
3070 double * opd_met_telfc_mcorr = NULL;
3071 double * opd_met_fc_corr = NULL;
3076 if ( !cpl_table_has_column (vis_SC,
"OPD_DISP") ) {
3077 cpl_msg_info (cpl_func,
"Cannot compute IMAGING_REF, not column OPD_DISP");
3079 return CPL_ERROR_NONE;
3082 cpl_array ** opd_disp = cpl_table_get_data_array (vis_SC,
"OPD_DISP");
3085 cpl_array ** phase_met_telfc = cpl_table_get_data_array (vis_SC,
"PHASE_MET_TELFC");
3094 cpl_array ** imaging_ref = cpl_table_get_data_array (vis_SC,
"IMAGING_REF");
3099 cpl_msg_info (cpl_func,
"imaging-ref-met = %s", imaging_ref_met);
3100 if (!strcmp (imaging_ref_met,
"TEL")) {
3101 cpl_msg_info (cpl_func,
"Use telescope metrology for IMAGING_REF computation");
3102 opd_met_telfc_mcorr = cpl_table_get_data_double (vis_SC,
"OPD_MET_TELFC_MCORR");
3103 opd_met_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
3104 }
else if (!strcmp (imaging_ref_met,
"FC_CORR")) {
3105 cpl_msg_info (cpl_func,
"Use corrected fiber coupler metrology for IMAGING_REF computation");
3106 opd_met_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
3107 }
else if (!strcmp (imaging_ref_met,
"FC")) {
3108 cpl_msg_info (cpl_func,
"Use fiber coupler metrology for IMAGING_REF computation");
3110 cpl_msg_error (cpl_func,
"Unknown metrology source for IMAGING_REF calculation!");
3111 cpl_ensure_code (imaging_ref_met, CPL_ERROR_ILLEGAL_INPUT);
3115 for (cpl_size row = 0; row < nrow; row ++) {
3118 imaging_ref[row] = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
3122 double opd_met_corr = (opd_met_telfc_mcorr ? opd_met_telfc_mcorr[row] : 0.0)
3123 + (opd_met_fc_corr ? opd_met_fc_corr[row] : 0.0);
3127 for (
int w = 0; w < nwave; w++) {
3129 double wavelength = cpl_table_get (wave_table,
"EFF_WAVE", w, NULL);
3133 cpl_array_set (imaging_ref[row], w,
3134 cpl_array_get (phase_ref[row], w, NULL)
3135 - cpl_array_get (opd_disp[row], w, NULL) * CPL_MATH_2PI / wavelength
3136 + (ucoord[row] * sep_U + vcoord[row] * sep_V) * CPL_MATH_2PI / wavelength
3137 - opd_met_corr * CPL_MATH_2PI / wavelength);
3140 if (!strcmp (imaging_ref_met,
"TEL"))
3141 cpl_array_set (imaging_ref[row], w,
3142 cpl_array_get (phase_ref[row], w, NULL)
3143 - cpl_array_get (opd_disp[row], w, NULL) * CPL_MATH_2PI / wavelength
3144 + (ucoord[row] * sep_U + vcoord[row] * sep_V) * CPL_MATH_2PI / wavelength
3145 - cpl_array_get (phase_met_telfc[row], w, NULL) );
3155 return CPL_ERROR_NONE;
3176 const cpl_parameterlist * parlist)
3179 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
3181 int nbase = 6,
ntel = 4, nclo = 4;
3191 cpl_size window_width = cpl_parameter_get_int(cpl_parameterlist_find_const(
3192 parlist,
"gravity.vis.pfactor-window-length"));
3195 if ( npol_sc != npol_ft ) {
3196 gravi_msg_warning (
"FIXME",
"Not sure this function works with npol_sc != npol_ft");
3208 for (
int pol = 0; pol < CPL_MAX(npol_sc,npol_ft); pol++) {
3211 cpl_msg_info (cpl_func,
"Start polarisation %d over %d",pol+1, CPL_MAX(npol_sc,npol_ft));
3212 cpl_msg_info(cpl_func,
"Insname FT : %s, pol %d npol %d",
3214 cpl_msg_info(cpl_func,
"Insname SC : %s, pol %d npol %d",
3231 CPLCHECK_MSG (
"Cannot get the VIS_MET and FDDL tables");
3268 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
3269 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
3270 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
3273 double * total_flux_ft = cpl_table_get_data_double (flux_FT,
"TOTALFLUX");
3278 for (cpl_size row = 0; row < nrow_ft *
ntel; row ++) {
3279 total_flux_ft[row] = cpl_array_get_mean (flux_ft[row]) * nwave_ft;
3283 cpl_table_duplicate_column (flux_FT,
"TOTALFLUX_UNSMOOTHED", flux_FT,
"TOTALFLUX");
3317 vis_FT, oi_wavelengthft);
3325 for (
int base = 0; base < nbase; base++) {
3328 sprintf (qc_name,
"ESO QC VFACTOR%s_P%d AVG",
GRAVI_BASE_NAME[base], pol+1);
3331 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean v-factor");
3333 sprintf (qc_name,
"ESO QC PFACTOR%s_P%d AVG",
GRAVI_BASE_NAME[base], pol+1);
3336 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean p-factor");
3339 for (
int closure = 0; closure < nclo; closure++) {
3342 sprintf (qc_name,
"ESO QC P3FACTOR%s_P%d AVG",
GRAVI_CLO_NAME[closure], pol+1);
3345 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean p3-factor");
3382 "GDELAY", oi_wavelengthft);
3386 "GDELAY", oi_wavelengthsc);
3390 "GDELAY_FT", oi_wavelengthft);
3408 return CPL_ERROR_NONE;
3421 cpl_ensure_code(p2vmred_data, CPL_ERROR_NULL_INPUT);
3422 cpl_ensure_code(
plist, CPL_ERROR_NULL_INPUT);
3427 int nbase = 6, nclo = 4;
3432 cpl_size nrow = cpl_table_get_nrow (vis_SC) / nbase;
3434 for (
int pol = 0; pol < CPL_MAX(npol_sc,npol_ft); pol++) {
3435 for (
int base = 0; base < nbase; base++) {
3438 sprintf (qc_name,
"ESO QC VFACTOR%s_P%d AVG",
GRAVI_BASE_NAME[base], pol+1);
3439 const cpl_property *qc_val = cpl_propertylist_get_property_const(p2vmred_header, qc_name);
3440 cpl_propertylist_append_property(
plist, qc_val);
3442 sprintf (qc_name,
"ESO QC PFACTOR%s_P%d AVG",
GRAVI_BASE_NAME[base], pol+1);
3443 qc_val = cpl_propertylist_get_property_const(p2vmred_header, qc_name);
3444 cpl_propertylist_append_property(
plist, qc_val);
3448 for (
int closure = 0; closure < nclo; closure++) {
3451 sprintf (qc_name,
"ESO QC P3FACTOR%s_P%d AVG",
GRAVI_CLO_NAME[closure], pol+1);
3452 const cpl_property *qc_val = cpl_propertylist_get_property_const(p2vmred_header, qc_name);
3453 cpl_propertylist_append_property(
plist, qc_val);
3458 cpl_propertylist_append_long_long(
plist,
"NROW", nrow);
3459 return CPL_ERROR_NONE;
3485 const cpl_parameterlist * parlist)
3488 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
3498 cpl_array * rejected_array = cpl_array_new(nbase,CPL_TYPE_DOUBLE);
3499 cpl_array_fill_window_double (rejected_array, 0, nbase, 0.0);
3505 cpl_msg_info (cpl_func,
"Cannot compute rejection flags for FT (no FT data)");
3508 cpl_msg_info (cpl_func,
"Compute rejection flags for FT");
3516 cpl_msg_info (cpl_func,
"SNR threshold to define fringe-detection in FT: %g", threshold_SNR_ft);
3517 cpl_msg_info (cpl_func,
"STATE threshold for FT: %g", threshold_STATE_ft);
3518 cpl_msg_info (cpl_func,
"GLOBAL_STATE threshold for FT: %g - %g", min_GSTATE_ft, max_GSTATE_ft);
3522 cpl_size nrow_ft = 0;
3523 for (
int pol = 0; pol < npol_ft; pol++) {
3527 double * snr = cpl_table_get_data_double (vis_FT,
"SNR_BOOT");
3528 int * state = cpl_table_get_data_int (vis_FT,
"STATE");
3529 int * gstate = cpl_table_get_data_int (vis_FT,
"OPDC_STATE");
3532 int * reject_flag_ft = cpl_table_get_data_int (vis_FT,
"REJECTION_FLAG");
3537 nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
3538 for (cpl_size row = 0; row < nrow_ft * nbase; row ++) {
3543 if ( snr[row] < threshold_SNR_ft )
3555 if ( state[row] < threshold_STATE_ft ||
3556 gstate[row] < min_GSTATE_ft ||
3557 gstate[row] > max_GSTATE_ft )
3568 cpl_array_set_double(rejected_array,row%nbase,
3569 cpl_array_get_double(rejected_array,row%nbase,NULL)+reject );
3575 cpl_array_multiply_scalar (rejected_array,100./(npol_ft*nrow_ft));
3578 for (
int base = 0; base < nbase; base++) {
3580 sprintf (qc_name,
"ESO QC REJECTED RATIO FT%s",
GRAVI_BASE_NAME[base]);
3581 double ratio = cpl_array_get_double (rejected_array, base, NULL);
3583 cpl_propertylist_update_int (p2vmred_header, qc_name, round(ratio));
3584 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[%] ratio of FT data flagged");
3587 sprintf (qc_name,
"ESO QC REJECTED RATIO FT");
3588 double ratio = cpl_array_get_mean (rejected_array);
3590 cpl_propertylist_update_int (p2vmred_header, qc_name, round(ratio));
3591 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[%] ratio of FT data flagged");
3592 cpl_array_fill_window_double (rejected_array, 0, nbase, 0.0);
3601 cpl_msg_info (cpl_func,
"Cannot compute tracking ratio for SC (no FT or SC data)");
3604 cpl_msg_info (cpl_func,
"Compute tracking ratio for each SC DIT");
3609 for (
int pol = 0; pol < npol_sc; pol++) {
3622 cpl_msg_info (cpl_func,
"Cannot compute rejection flags for SC (no SC data)");
3625 cpl_msg_info (cpl_func,
"Compute rejection flags for SC");
3631 cpl_msg_info (cpl_func,
"Fringe-detection ratio to discard frame on SC: %g", minlockratio_sc);
3632 cpl_msg_info (cpl_func,
"vFactor threshold to discard frame on SC: %g", minvfactor_sc);
3635 double opd_pupil_stddev_max_sc =
gravi_param_get_double (parlist,
"gravity.signal.opd-pupil-stddev-max-sc");
3637 cpl_msg_info (cpl_func,
"OPD_PUPIL threshold (abs) to discard frame on SC: %g", opd_pupil_max_sc);
3638 cpl_msg_info (cpl_func,
"OPD_PUPIL_STDDEV threshold to discard frame on SC: %g", opd_pupil_stddev_max_sc);
3642 cpl_size nrow_sc = 0;
3643 for (
int pol = 0; pol < npol_sc; pol++) {
3647 double * vFactor_wl = cpl_table_get_data_double (vis_SC,
"V_FACTOR_WL");
3648 double * fringedet_ftdit = cpl_table_get_data_double (vis_SC,
"FRINGEDET_RATIO");
3652 double * opd_metdit_pupil = NULL;
3653 if ( cpl_table_has_column (vis_SC,
"OPD_MET_PUPIL") )
3654 opd_metdit_pupil = cpl_table_get_data_double (vis_SC,
"OPD_MET_PUPIL");
3656 double * opd_metdit_pupil_stddev = NULL;
3657 if ( cpl_table_has_column (vis_SC,
"OPD_MET_PUPIL_STDDEV") )
3658 opd_metdit_pupil_stddev = cpl_table_get_data_double (vis_SC,
"OPD_MET_PUPIL_STDDEV");
3663 int * reject_flag_sc = cpl_table_get_data_int (vis_SC,
"REJECTION_FLAG");
3668 nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
3669 for (cpl_size row_sc = 0; row_sc < nrow_sc * nbase; row_sc++) {
3674 if ( fringedet_ftdit[row_sc] < minlockratio_sc )
3683 int vfactor_bit = 1;
3684 if ( vFactor_wl[row_sc] < minvfactor_sc )
3695 int opd_pupil_bit = 2;
3696 if ( opd_metdit_pupil ) {
3697 if ( fabs(opd_metdit_pupil[row_sc]) > opd_pupil_max_sc )
3709 int opd_pupil_stddev_bit = 3;
3710 if ( opd_metdit_pupil_stddev ) {
3711 if ( opd_metdit_pupil_stddev[row_sc] > opd_pupil_stddev_max_sc )
3713 gravi_bit_set (reject_flag_sc[row_sc], opd_pupil_stddev_bit);
3723 cpl_array_set_double(rejected_array,row_sc%nbase,
3724 cpl_array_get_double(rejected_array,row_sc%nbase,NULL)+reject );
3729 cpl_array_multiply_scalar (rejected_array,100./(npol_sc*nrow_sc));
3732 for (
int base = 0; base < nbase; base++) {
3734 sprintf (qc_name,
"ESO QC REJECTED RATIO SC%s",
GRAVI_BASE_NAME[base]);
3735 double ratio = cpl_array_get_double (rejected_array, base, NULL);
3737 cpl_propertylist_update_int (p2vmred_header, qc_name, ratio);
3738 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[%] ratio of SC data flagged");
3741 sprintf (qc_name,
"ESO QC REJECTED RATIO SC");
3742 double ratio = cpl_array_get_mean (rejected_array);
3744 cpl_propertylist_update_int (p2vmred_header, qc_name, ratio);
3745 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[%] ratio of SC data flagged");
3748 cpl_array_delete(rejected_array);
3750 return CPL_ERROR_NONE;
#define gravi_table_get_value(table, name, row, value)
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
#define gravi_data_get_oi_t3(data, type, pol, npol)
#define gravi_data_get_oi_flux(data, type, pol, npol)
#define gravi_data_get_header(data)
#define gravi_data_get_oi_wave(data, type, pol, npol)
#define gravi_data_get_oi_vis(data, type, pol, npol)
cpl_msg_debug(cpl_func, "Spectra has <50 pixels -> don't flat")
cpl_propertylist * header
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
cpl_propertylist_update_double(header, "ESO QC MINWAVE SC", cpl_propertylist_get_double(plist, "ESO QC MINWAVE SC"))
#define GRAVI_OI_VIS_ACQ_EXT
#define GRAVI_INSNAME(type, pol, npol)
#define GRAVI_OI_VIS_MET_EXT
#define CPLCHECK_INT(msg)
#define gravi_msg_function_exit(flag)
#define gravi_bit_clear(number, pos)
#define FREE(function, variable)
#define gravi_msg_function_start(flag)
#define CPLCHECK_MSG(msg)
#define gravi_bit_set(number, pos)
#define FREELOOP(function, variable, n)
cpl_error_code gravi_table_smooth_column(cpl_table *oi_vis, const char *input_name, const char *output_name, int nsmooth, int nbase)
cpl_error_code gravi_table_add_columns(cpl_table *oi_vis1, const char *name1, cpl_table *oi_vis2, const char *name2)
cpl_error_code gravi_table_runint_column(cpl_table *oi_vis, const char *input_name, const char *output_name, int nsmooth, int nbase)
double gravi_table_get_column_mean(cpl_table *table, const char *name, int base, int nbase)
cpl_error_code gravi_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
cpl_error_code gravi_table_compute_group_delay(cpl_table *table, const char *input, const char *flag, const char *output, cpl_table *oi_wave)
cpl_error_code gravi_table_init_column_array(cpl_table *table, const char *name, const char *unit, cpl_type type, cpl_size size)
cpl_error_code gravi_array_add_phase(cpl_array *input, double factor, cpl_array *phase)
Add a REAL phase to a REAL phase array, in-place: input = input + factor * phase.
cpl_error_code gravi_vector_abs(cpl_vector *vector)
Return the running median of a vector, with special care for the boundaray, that are filled with the ...
cpl_array * gravi_array_init_double(long n, double value)
double ** gravi_table_get_data_array_double(cpl_table *table, const char *name)
cpl_vector * gravi_vector_median(const cpl_vector *vector, cpl_size hw)
Return the running median of a vector, with special care for the boundaray, that are filled with the ...
cpl_array * gravi_table_get_column_sum_array(cpl_table *table, const char *name, int base, int nbase)
cpl_error_code gravi_array_phase_wrap(cpl_array *input)
cpl_array * gravi_array_init_double_complex(long n, double complex value)
cpl_error_code gravi_array_multiply_phasor(cpl_array *input, double complex factor, cpl_array *phase)
Multiply a REAL phase to a COMPLEX array, in-place: input = input * cexp (factor * phase)
cpl_error_code gravi_table_new_column_array(cpl_table *table, const char *name, const char *unit, cpl_type type, cpl_size size)
int gravi_data_has_extension(gravi_data *raw_calib, const char *ext_name)
Check if data has extension with given EXTNAME.
int gravi_data_has_type(gravi_data *self, const char *type)
Return the number of ext whose EXTNAME and INSNAME match 'type'.
cpl_table * gravi_data_get_table(gravi_data *self, const char *extname)
Return a pointer on a table extension by its EXTNAME.
int gravi_param_get_bool(const cpl_parameterlist *parlist, const char *name)
const char * gravi_param_get_string(const cpl_parameterlist *parlist, const char *name)
int gravi_param_get_int(const cpl_parameterlist *parlist, const char *name)
double gravi_param_get_double(const cpl_parameterlist *parlist, const char *name)
int gravi_pfits_get_pola_num(const cpl_propertylist *plist, int type_data)
double gravi_pfits_get_gdzero(const cpl_propertylist *plist, int tel)
double gravi_pfits_get_sobj_y(const cpl_propertylist *plist)
double gravi_pfits_get_oplzero(const cpl_propertylist *plist, int tel)
double gravi_pfits_get_dit_sc(const cpl_propertylist *plist)
int gravi_pfits_has_gdzero(const cpl_propertylist *plist, int tel)
double gravi_pfits_get_period_sc(const cpl_propertylist *plist)
double gravi_pfits_get_sobj_x(const cpl_propertylist *plist)
cpl_error_code gravi_vis_create_f1f2_sc(cpl_table *vis_SC, cpl_table *flux_SC)
Compute the photometric normalisation for the SC.
cpl_error_code gravi_vis_create_pfactor_sc(cpl_table *vis_SC, cpl_table *t3_SC, cpl_table *flux_FT)
Compute the PFACTOR for the SC.
cpl_error_code gravi_vis_create_vfactor_sc(cpl_table *vis_SC, cpl_table *wave_table_sc, cpl_table *vis_FT, cpl_table *wave_table_ft)
Compute the VFACTOR for each SC DIT based on real-time FT.
cpl_error_code gravi_flux_create_acq_sc(cpl_table *vis_SC, cpl_table *vis_ACQ)
Compute the averaged ACQ signal for each SC DIT per beam.
cpl_error_code gravi_vis_compute_isdelay(cpl_table *oi_vis, const char *name_isp, const char *name_gdl, cpl_table *oi_wavelength)
Compute the group-delay from interspectra.
cpl_error_code gravi_vis_create_f1f2_ft(cpl_table *vis_FT, cpl_table *flux_FT)
Compute the photometric normalisation for the FT.
cpl_error_code gravi_vis_create_lockratio_sc(cpl_table *vis_SC, cpl_table *vis_FT)
Create the clockratio for each SC DIT and baseline.
cpl_error_code gravi_compute_snr(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Compute real-time SNR and Group-Delay of the observation.
cpl_error_code gravi_flux_create_fddlpos_sc(cpl_table *flux_SC, cpl_table *fddl_table)
Compute the averaged FDDL signal for each SC DIT per beam.
cpl_error_code gravi_flux_create_fddllin_sc(cpl_table *flux_SC, cpl_table *disp_table)
Compute the (FDDL_SC + FDDL_FT)/2 position in [m].
cpl_error_code gravi_vis_create_opdsc_ft(cpl_table *vis_FT, cpl_table *vis_SC, double dit_sc)
Compute the resampled SC signal for each FT DIT per base.
cpl_error_code gravi_flux_create_totalflux_sc(cpl_table *flux_SC, cpl_table *flux_FT)
Compute total flux of each DIT for the SC and of the FT.
cpl_error_code gravi_flux_create_met_sc(cpl_table *flux_SC, cpl_table *vis_MET)
Compute the averaged MET signal for each SC DIT per beam.
cpl_error_code gravi_vis_correct_phasediff(cpl_table *oi_vis1, const char *name1, cpl_table *oi_vis2, const char *name2, double *phasediff)
Correct for mean phase-difference between coherent fluxes.
cpl_error_code gravi_vis_create_phaseref_ft(cpl_table *vis_FT)
Compute the self-reference phase for each FT DIT.
cpl_error_code gravi_vis_create_pfactor_ft(cpl_table *vis_FT, cpl_table *flux_FT, cpl_size window_width)
Compute the PFACTOR for the FT.
cpl_error_code gravi_vis_create_acq_sc(cpl_table *vis_SC, cpl_table *vis_ACQ)
Compute the averaged ACQ signal for each SC DIT per base.
cpl_error_code gravi_vis_create_opddisp_sc(cpl_table *vis_SC, cpl_table *flux_SC, cpl_table *wave_table, cpl_table *disp_table, cpl_propertylist *header, const cpl_parameterlist *parlist)
Compute the MET opd including the dispersion from the FDDL.
cpl_error_code gravi_vis_compute_interspectre(cpl_table *oi_vis, const char *name_vis, const char *name_is, const char *name_flag)
Compute the real-time interspectra.
cpl_error_code gravi_vis_create_met_ft(cpl_table *vis_FT, cpl_table *vis_MET)
Compute the resampled MET signal for each FT DIT per base.
cpl_error_code gravi_vis_create_phaseref_sc(cpl_table *vis_SC, cpl_table *wave_table_sc, cpl_table *wave_table_ft, cpl_propertylist *header, const cpl_parameterlist *parlist)
Compute the reference phase for each SC DIT.
cpl_error_code gravi_compute_signals(gravi_data *p2vmred_data, gravi_data *disp_data, const cpl_parameterlist *parlist)
Create intermediate signal in the P2VMREDUCED file.
cpl_error_code gravi_vis_create_imagingref_sc(cpl_table *vis_SC, cpl_table *wave_table, cpl_propertylist *header, const cpl_parameterlist *parlist)
Create phase-referenced imaging data in the P2VMREDUCED file.
cpl_error_code gravi_signal_create_sync(cpl_table *vis_SC, int nbase_sc, double dit_sc, cpl_table *vis_FT, int nbase_ft, const char *name)
Compute synchronisation indices between OIFITS tables.
cpl_error_code gravi_vis_create_met_sc(cpl_table *vis_SC, cpl_table *vis_MET, cpl_table *wave_table)
Compute the averaged MET signal for each SC DIT per base.
cpl_error_code gravi_compute_outliers(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Compute the outliers flags.
cpl_error_code gravi_compute_rejection(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Create rejection flags P2VMREDUCED file.
cpl_error_code gravi_vis_compute_snr(cpl_table *oi_vis, const char *name_pha, const char *name_var, const char *name_snr)
Compute the real-time SNR.
cpl_error_code gravi_create_outlier_flag_ft(cpl_table *flux_SC, cpl_table *vis_SC)
Create the list of outlier. For the FT, this is filled with 0.
cpl_error_code gravi_copy_p2vm_qcs(gravi_data *p2vmred_data, cpl_propertylist *plist)
Copy PFACTOR and VFACTOR QCs so that they may be aggregated over all frames.
cpl_error_code gravi_vis_bootstrap_snr_and_delay(cpl_table *oi_vis, const char *name_snr, const char *name_gdl)
Boostrap real-time SNR and GDELAY.
cpl_error_code gravi_vis_compute_mean_phasor(cpl_table *oi_vis, const char *name_vis, const char *name_err, const char *name_pha, const char *name_var, const char *name_flag)
Compute real-time mean phasor of a VISDATA by averaging all spectral elements.
cpl_error_code gravi_create_outlier_flag_sc(cpl_table *flux_SC, cpl_table *vis_SC, double chi2r_threshold, double chi2r_sigma)
Create the list of outlier based on the values of the chi2.
int GRAVI_TRI_SIGN[GRAVI_NBASE][2][2]
cpl_error_code gravi_msg_warning(const char *component, const char *msg)
char GRAVI_BASE_NAME[6][3]
int GRAVI_BASE_TEL[GRAVI_NBASE][2]
char GRAVI_CLO_NAME[4][4]
int GRAVI_TRI_BASE[GRAVI_NBASE][2][2]