52 #include "gravi_data.h" 53 #include "gravi_dfs.h" 54 #include "gravi_pfits.h" 55 #include "gravi_cpl.h" 57 #include "gravi_utils.h" 58 #include "gravi_signal.h" 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);
76 const char * name_vis,
77 const char * name_is);
79 const char * name_pha,
80 const char * name_var,
81 const char * name_snr);
83 const char * name_isp,
84 const char * name_gdl,
85 cpl_table * oi_wavelength);
101 cpl_table * wave_table_sc,
103 cpl_table * wave_table_ft);
108 cpl_table * wave_table_sc,
109 cpl_table * wave_table_ft);
114 cpl_table * wave_table,
115 cpl_table * disp_table,
116 cpl_propertylist * header,
117 const cpl_parameterlist * parlist);
120 cpl_table * wave_table,
121 cpl_propertylist * header,
122 const cpl_parameterlist * parlist);
144 const char * name_snr,
145 const char * name_gdl)
147 gravi_msg_function_start(0);
148 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
149 cpl_ensure_code (name_snr, CPL_ERROR_NULL_INPUT);
150 cpl_ensure_code (name_gdl, CPL_ERROR_NULL_INPUT);
153 double * snr = cpl_table_get_data_double (oi_vis, name_snr);
154 double * gdl = cpl_table_get_data_double (oi_vis, name_gdl);
156 CPLCHECK_MSG(
"Cannot load data to compute SNR column");
159 cpl_size nrow = cpl_table_get_nrow (oi_vis) / 6;
160 for (
int base=0; base<6; base++ ) {
163 for (
int tri=0; tri<2; tri++ ) {
166 int b1 = GRAVI_TRI_BASE[base][tri][0];
167 int b2 = GRAVI_TRI_BASE[base][tri][1];
168 int sign1 = GRAVI_TRI_SIGN[base][tri][0];
169 int sign2 = GRAVI_TRI_SIGN[base][tri][1];
171 cpl_msg_debug(cpl_func,
"Found triangle tels %i%i -> bases (%i,%i, %i,%i)",
172 GRAVI_BASE_TEL[base][0],GRAVI_BASE_TEL[base][1],b1,b2,sign1,sign2);
175 for (cpl_size row = 0; row < nrow; row++) {
177 double snrN = CPL_MIN( snr[row*6+b1], snr[row*6+b2] );
178 if ( snrN > snr[row*6+base] ) {
179 snr[row*6+base] = snrN;
180 gdl[row*6+base] = sign1 * gdl[row*6+b1] + sign2 * gdl[row*6+b2];
187 cpl_msg_debug(cpl_func,
"TODO: loop on 3-baseline bootstrap");
191 gravi_msg_function_exit(0);
192 return CPL_ERROR_NONE;
212 cpl_table * oi_vis2,
const char *name2,
215 gravi_msg_function_start(0);
216 cpl_ensure_code (oi_vis1, CPL_ERROR_NULL_INPUT);
217 cpl_ensure_code (oi_vis2, CPL_ERROR_NULL_INPUT);
218 cpl_ensure_code (name1, CPL_ERROR_NULL_INPUT);
219 cpl_ensure_code (name2, CPL_ERROR_NULL_INPUT);
220 cpl_ensure_code (phasediff, CPL_ERROR_ILLEGAL_OUTPUT);
224 cpl_type type1 = cpl_table_get_column_type (oi_vis1, name1 );
225 cpl_type type2 = cpl_table_get_column_type (oi_vis2, name2 );
226 cpl_size nrow1 = cpl_table_get_nrow (oi_vis1) / nbase;
227 cpl_size nrow2 = cpl_table_get_nrow (oi_vis2) / nbase;
229 if ( type1 != CPL_TYPE_DOUBLE_COMPLEX ||
230 type2 != CPL_TYPE_DOUBLE_COMPLEX ||
232 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
"input columns not conformables or not DOUBLE COMPLEX");
234 double complex * data1 = cpl_table_get_data_double_complex (oi_vis1, name1);
235 double complex * data2 = cpl_table_get_data_double_complex (oi_vis2, name2);
236 CPLCHECK_MSG(
"Cannot load data");
239 for (
int base = 0; base < nbase; base++) {
240 double complex phasor = 0.0;
243 for (cpl_size row = 0; row < nrow1; row++) {
244 phasor += data1[row*nbase+base] * conj( data2[row*6+base] );
247 phasediff[base] = carg( phasor );
250 for (cpl_size row = 0; row < nrow1; row++) {
251 data1[row*6+base] *= cexp( -I*phasediff[base] );
255 gravi_msg_function_exit(0);
256 return CPL_ERROR_NONE;
273 const char * name_vis,
274 const char * name_err,
275 const char * name_pha,
276 const char * name_var)
278 gravi_msg_function_start(0);
279 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
280 cpl_ensure_code (name_vis, CPL_ERROR_NULL_INPUT);
281 cpl_ensure_code (name_err, CPL_ERROR_NULL_INPUT);
282 cpl_ensure_code (name_pha, CPL_ERROR_ILLEGAL_OUTPUT);
283 cpl_ensure_code (name_var, CPL_ERROR_ILLEGAL_OUTPUT);
285 cpl_size nrow = cpl_table_get_nrow (oi_vis);
288 gravi_table_new_column (oi_vis, name_pha,
"e", CPL_TYPE_DOUBLE_COMPLEX);
289 double complex * phasorRaw = cpl_table_get_data_double_complex (oi_vis, name_pha);
291 gravi_table_new_column (oi_vis, name_var,
"e^2", CPL_TYPE_DOUBLE);
292 double * varRaw = cpl_table_get_data_double (oi_vis, name_var);
295 cpl_array** tVis = cpl_table_get_data_array (oi_vis, name_vis);
296 cpl_array** tErr = cpl_table_get_data_array (oi_vis, name_err);
299 int nwave = cpl_array_get_size (tVis[0]);
301 CPLCHECK_MSG (
"Cannot get data");
304 for (cpl_size row = 0; row < nrow; row++ ) {
306 double complex * cpx = cpl_array_get_data_double_complex (tVis[row]);
307 double complex * err = cpl_array_get_data_double_complex (tErr[row]);
308 cpl_ensure_code (cpx, CPL_ERROR_ILLEGAL_INPUT);
309 cpl_ensure_code (err, CPL_ERROR_ILLEGAL_INPUT);
310 phasorRaw[row] = 0.0;
314 for (cpl_size wave = 0; wave < nwave; wave++) {
315 phasorRaw[row] += cpx[wave];
316 varRaw[row] += cabs(err[wave]) * cabs(err[wave]);
322 CPLCHECK_MSG(
"Cannot fill IS and PHASOR column");
324 gravi_msg_function_exit(0);
325 return CPL_ERROR_NONE;
342 const char * name_vis,
343 const char * name_is)
345 gravi_msg_function_start(0);
346 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
347 cpl_ensure_code (name_vis, CPL_ERROR_NULL_INPUT);
348 cpl_ensure_code (name_is, CPL_ERROR_NULL_INPUT);
350 cpl_size nrow = cpl_table_get_nrow (oi_vis);
353 cpl_array** tVis = cpl_table_get_data_array (oi_vis, name_vis);
354 cpl_ensure_code (tVis, CPL_ERROR_ILLEGAL_INPUT);
357 int nwave = cpl_array_get_size (tVis[0]);
359 CPLCHECK_MSG (
"Cannot get data");
362 gravi_table_new_column (oi_vis, name_is,
"e^2", CPL_TYPE_DOUBLE_COMPLEX);
363 double complex * interSpectreRaw = cpl_table_get_data_double_complex (oi_vis, name_is);
365 CPLCHECK_MSG (
"Cannot create columns");
368 for (cpl_size row = 0; row < nrow; row++) {
370 double complex * cpx = cpl_array_get_data_double_complex (tVis[row]);
371 cpl_ensure_code (cpx, CPL_ERROR_ILLEGAL_INPUT);
373 interSpectreRaw[row] = 0.0;
376 for (cpl_size wave = 0; wave < nwave-1; wave++)
377 interSpectreRaw[row] += cpx[wave] * conj( cpx[wave+1] );
381 CPLCHECK_MSG (
"Cannot fill IS column");
383 gravi_msg_function_exit(0);
384 return CPL_ERROR_NONE;
402 const char * name_pha,
403 const char * name_var,
404 const char * name_snr)
406 gravi_msg_function_start(0);
407 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
408 cpl_ensure_code (name_pha, CPL_ERROR_NULL_INPUT);
409 cpl_ensure_code (name_var, CPL_ERROR_NULL_INPUT);
410 cpl_ensure_code (name_snr, CPL_ERROR_ILLEGAL_OUTPUT);
412 cpl_size nrow = cpl_table_get_nrow (oi_vis);
415 double complex * cpx = cpl_table_get_data_double_complex (oi_vis, name_pha);
416 double * dbl = cpl_table_get_data_double (oi_vis, name_var);
418 CPLCHECK_MSG (
"Cannot load data to compute SNR column");
421 gravi_table_new_column (oi_vis, name_snr, NULL, CPL_TYPE_DOUBLE);
422 double * snr = cpl_table_get_data_double (oi_vis, name_snr);
424 CPLCHECK_MSG (
"Cannot createl SNR column");
427 for (cpl_size row = 0; row < nrow; row++) {
428 snr[row] = cabs( cpx[row] ) / sqrt (fabs(dbl[row]));
431 CPLCHECK_MSG (
"Cannot fill SNR column");
433 gravi_msg_function_exit(0);
434 return CPL_ERROR_NONE;
452 const char * name_isp,
453 const char * name_gdl,
454 cpl_table * oi_wavelength)
456 gravi_msg_function_start(0);
457 cpl_ensure_code (oi_vis, CPL_ERROR_NULL_INPUT);
458 cpl_ensure_code (name_isp, CPL_ERROR_NULL_INPUT);
459 cpl_ensure_code (name_gdl, CPL_ERROR_NULL_INPUT);
460 cpl_ensure_code (oi_wavelength, CPL_ERROR_NULL_INPUT);
462 cpl_size nrow = cpl_table_get_nrow (oi_vis);
465 double complex * cpx1 = cpl_table_get_data_double_complex (oi_vis, name_isp);
467 CPLCHECK_MSG(
"Cannot load data to compute GDELAY column");
470 gravi_table_new_column (oi_vis, name_gdl,
"m", CPL_TYPE_DOUBLE);
471 double * gdl = cpl_table_get_data_double (oi_vis, name_gdl);
472 CPLCHECK_MSG(
"Cannot set GDELAY column");
475 cpl_size wave = cpl_table_get_ncol (oi_wavelength)/2;
476 double factor = 1./fabs(1./cpl_table_get (oi_wavelength,
"EFF_WAVE", wave, NULL) -
477 1./cpl_table_get (oi_wavelength,
"EFF_WAVE", wave+1, NULL));
479 cpl_msg_debug (cpl_func,
"Compute the coherence length (%e m)", factor);
481 CPLCHECK_MSG(
"Cannot compute coherence length");
484 for (cpl_size row = 0; row < nrow; row++ ) {
486 gdl[row] = (double)carg( cpx1[row] ) / (2.0*M_PI) * factor;
489 CPLCHECK_MSG(
"Cannot fill GDELAY column");
491 gravi_msg_function_exit(0);
492 return CPL_ERROR_NONE;
512 const cpl_parameterlist * parlist)
515 gravi_msg_function_start(1);
516 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
521 cpl_propertylist * p2vmred_header = gravi_data_get_header (p2vmred_data);
522 double periode_sc = gravi_pfits_get_period_sc (p2vmred_header);
523 CPLCHECK_MSG (
"Cannot get header");
531 for (
int type_data = 0; type_data < 2; type_data ++) {
534 cpl_msg_info (cpl_func,
"Cannot compute real-time snr for %s " 535 "(no data)", GRAVI_TYPE(type_data));
540 if (type_data == GRAVI_FT) {
543 nsmooth = gravi_param_get_int (parlist,
"gravity.signal.nsmooth-snr-ft");
548 nsmooth = 1./ periode_sc;
549 nsmooth = CPL_MIN (CPL_MAX (nsmooth, 0), 20);
554 int npol = gravi_pfits_get_pola_num (p2vmred_header, type_data);
555 for (
int pol = 0; pol < npol; pol++) {
557 cpl_msg_info(cpl_func,
"Compute SNR for pol %i of %s (smoothed over %i frames)",
558 pol+1,GRAVI_TYPE(type_data),2*nsmooth+1);
560 cpl_table * oi_vis = gravi_data_get_oi_vis (p2vmred_data, type_data, pol, npol);
561 CPLCHECK_MSG (
"Cannot get data");
569 CPLCHECK_MSG (
"Cannot compute PHASOR and IPHASOR");
573 gravi_table_runint_column (oi_vis,
"IPHASOR",
"IPHASOR_SMT", nsmooth, nbase);
574 gravi_table_runint_column (oi_vis,
"PHASOR",
"PHASOR_SMT", nsmooth, nbase);
575 gravi_table_runint_column (oi_vis,
"PHASOR_VAR",
"PHASOR_VAR_SMT", nsmooth, nbase);
576 CPLCHECK_MSG (
"Cannot running integration");
588 for (
int type_data = 0; type_data < 2; type_data ++) {
591 cpl_msg_info (cpl_func,
"Cannot compute bootstraped snr %s " 592 "(no data)", GRAVI_TYPE(type_data));
596 int npol = gravi_pfits_get_pola_num (p2vmred_header, type_data);
598 cpl_msg_info(cpl_func,
"Compute bootstraped GDELAY_BOOT and SNR for %s...",(type_data==GRAVI_FT?
"FT":
"SC"));
601 cpl_table * oi_wave = gravi_data_get_oi_wave (p2vmred_data, type_data, 0, npol);
602 cpl_table * oi_vis_p2, * oi_vis_p1;
603 oi_vis_p1 = gravi_data_get_oi_vis (p2vmred_data, type_data, 0, npol);
605 CPLCHECK_MSG(
"Cannot get OI_VIS_P1 table");
608 cpl_table_duplicate_column (oi_vis_p1,
"IPHASOR_BOOT", oi_vis_p1,
"IPHASOR_SMT");
609 CPLCHECK_MSG(
"Cannot duplicate columns");
612 cpl_table_duplicate_column (oi_vis_p1,
"PHASOR_BOOT", oi_vis_p1,
"PHASOR_SMT");
613 CPLCHECK_MSG(
"Cannot duplicate columns");
615 cpl_table_duplicate_column (oi_vis_p1,
"PHASOR_VAR_BOOT", oi_vis_p1,
"PHASOR_VAR_SMT");
616 CPLCHECK_MSG(
"Cannot duplicate columns");
621 cpl_msg_info(cpl_func,
"Add the signal of both polarisation to enhance SNR and GDELAY accuracy");
622 oi_vis_p2 = gravi_data_get_oi_vis (p2vmred_data, type_data, 1, npol);
624 CPLCHECK_MSG(
"Cannot get OI_VIS_P2 table");
628 gravi_table_add_columns (oi_vis_p1,
"IPHASOR_BOOT", oi_vis_p2,
"IPHASOR_SMT");
629 CPLCHECK_MSG(
"Cannot add columns");
637 gravi_table_add_columns (oi_vis_p1,
"PHASOR_BOOT", oi_vis_p2,
"PHASOR_SMT");
638 gravi_table_add_columns (oi_vis_p1,
"PHASOR_VAR_BOOT", oi_vis_p2,
"PHASOR_VAR_SMT");
641 cpl_msg_info(cpl_func,
"Add the phase difference as QC parameters");
642 for (
int base = 0; base < nbase; base++) {
643 sprintf (qc_name,
"ESO QC PHASE_POLDIFF_%s%d%d",(type_data==GRAVI_FT?
"FT":
"SC"),
644 GRAVI_BASE_TEL[base][0]+1, GRAVI_BASE_TEL[base][1]+1);
645 cpl_propertylist_update_double (p2vmred_header, qc_name, phasediff[base] * CPL_MATH_DEG_RAD);
646 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[deg] differential phase" );
647 cpl_msg_info (cpl_func,
"%s=%f [deg]",qc_name,phasediff[base] * CPL_MATH_DEG_RAD);
648 CPLCHECK_MSG(
"QC PHASE_POLDIFF");
651 CPLCHECK_MSG(
"Cannot add columns of second polarisation");
661 int nstabilize = (type_data == GRAVI_FT ? 50 : 1);
662 gravi_table_smooth_column (oi_vis_p1,
"SNR_BOOT_TMP",
"SNR_BOOT", nstabilize, nbase);
665 CPLCHECK_MSG(
"Cannot compute and fill SNR_BOOT or GDELAY_BOOT column");
669 cpl_msg_info (cpl_func,
"Duplicate columns in polarisation 2");
670 if (cpl_table_has_column (oi_vis_p2,
"SNR_BOOT"))
671 cpl_table_erase_column(oi_vis_p2,
"SNR_BOOT");
672 cpl_table_duplicate_column (oi_vis_p2,
"SNR_BOOT", oi_vis_p1,
"SNR_BOOT");
674 if (cpl_table_has_column (oi_vis_p2,
"GDELAY_BOOT"))
675 cpl_table_erase_column(oi_vis_p2,
"GDELAY_BOOT");
676 cpl_table_duplicate_column (oi_vis_p2,
"GDELAY_BOOT", oi_vis_p1,
"GDELAY_BOOT");
677 CPLCHECK_MSG(
"Cannot duplicate column for polarisation 2");
682 cpl_msg_info (cpl_func,
"Erase useless columns in P2VMREDUCED");
684 cpl_table_erase_column (oi_vis_p1,
"IPHASOR");
685 cpl_table_erase_column (oi_vis_p1,
"IPHASOR_SMT");
686 cpl_table_erase_column (oi_vis_p1,
"IPHASOR_BOOT");
687 cpl_table_erase_column (oi_vis_p1,
"PHASOR");
688 cpl_table_erase_column (oi_vis_p1,
"PHASOR_SMT");
689 cpl_table_erase_column (oi_vis_p1,
"PHASOR_BOOT");
690 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR");
691 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR_SMT");
692 cpl_table_erase_column (oi_vis_p1,
"PHASOR_VAR_BOOT");
693 cpl_table_erase_column (oi_vis_p1,
"SNR_BOOT_TMP");
694 CPLCHECK_MSG(
"Cannot erase column for polarisation 1");
697 cpl_table_erase_column (oi_vis_p2,
"IPHASOR");
698 cpl_table_erase_column (oi_vis_p2,
"IPHASOR_SMT");
699 cpl_table_erase_column (oi_vis_p2,
"PHASOR");
700 cpl_table_erase_column (oi_vis_p2,
"PHASOR_SMT");
701 cpl_table_erase_column (oi_vis_p2,
"PHASOR_VAR");
702 cpl_table_erase_column (oi_vis_p2,
"PHASOR_VAR_SMT");
703 CPLCHECK_MSG(
"Cannot erase column for polarisation 2");
710 for (cpl_size base = 0; base < nbase; base++) {
711 snprintf (qc_name, 100,
"ESO QC SNRB_%s%s AVG", GRAVI_TYPE(type_data), GRAVI_BASE_NAME[base]);
712 qc_value = gravi_table_get_column_mean (oi_vis_p1,
"SNR_BOOT", base, nbase);
713 cpl_propertylist_update_double (p2vmred_header, qc_name, qc_value);
714 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean bootstrapped SNR");
719 gravi_msg_function_exit(1);
720 return CPL_ERROR_NONE;
744 cpl_table * vis_FT,
int nbase_ft,
747 gravi_msg_function_start(1);
748 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
749 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
752 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase_sc;
753 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase_ft;
755 CPLCHECK_MSG (
"Cannot get data");
759 sprintf (full_name,
"FIRST_%s", name);
760 gravi_table_new_column (vis_SC, full_name, NULL, CPL_TYPE_INT);
761 int * first_ft = cpl_table_get_data_int (vis_SC, full_name);
763 sprintf (full_name,
"LAST_%s", name);
764 gravi_table_new_column (vis_SC, full_name, NULL, CPL_TYPE_INT);
765 int * last_ft = cpl_table_get_data_int (vis_SC, full_name);
767 sprintf (full_name,
"NFRAME_%s", name);
768 gravi_table_new_column (vis_SC, full_name, NULL, CPL_TYPE_INT);
769 int * nframe_ft = cpl_table_get_data_int (vis_SC, full_name);
771 CPLCHECK_MSG (
"Cannot create columns");
774 for (cpl_size base_sc = 0; base_sc < nbase_sc; base_sc++) {
779 cpl_size base_ft = base_sc % nbase_ft;
782 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
785 while ( cpl_table_get (vis_FT,
"TIME", row_ft * nbase_ft + base_ft, NULL) <
786 cpl_table_get (vis_SC,
"TIME", row_sc * nbase_sc + base_sc, NULL) - dit_sc/2.) {
787 if (row_ft >= nrow_ft-1)
break;
790 first_ft[row_sc * nbase_sc + base_sc] = row_ft;
793 while ( cpl_table_get (vis_FT,
"TIME", row_ft * nbase_ft + base_ft, NULL) <
794 cpl_table_get (vis_SC,
"TIME", row_sc * nbase_sc + base_sc, NULL) + dit_sc/2.) {
795 if (row_ft >= nrow_ft-1)
break;
798 last_ft[row_sc * nbase_sc + base_sc] = row_ft;
801 if ( first_ft[row_sc * nbase_sc + base_sc] < 2 ||
802 last_ft[row_sc * nbase_sc + base_sc] > nrow_ft - 2 ) {
803 cpl_msg_warning (cpl_func,
"Not enough %s data to synchronise with DIT %lli over %lli", name, row_sc+1, nrow_sc);
804 first_ft[row_sc * nbase_sc + base_sc] = 0;
805 last_ft[row_sc * nbase_sc + base_sc] = 0;
809 nframe_ft[row_sc * nbase_sc + base_sc] = last_ft[row_sc * nbase_sc + base_sc] - first_ft[row_sc * nbase_sc + base_sc];
814 gravi_msg_function_exit(1);
815 return CPL_ERROR_NONE;
832 gravi_msg_function_start(1);
833 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
834 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
837 cpl_size nbase = 6, ntel = 4;
838 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
841 gravi_table_new_column (vis_SC,
"P_FACTOR", NULL, CPL_TYPE_DOUBLE);
842 double * pFactor = cpl_table_get_data_double (vis_SC,
"P_FACTOR");
845 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
846 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
849 double * flux = cpl_table_get_data_double (flux_FT,
"TOTALFLUX");
851 CPLCHECK_MSG (
"Cannot get pointer to data");
854 for (cpl_size base = 0; base < nbase; base++) {
855 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
856 cpl_size nsc = row_sc*nbase+base;
859 double sf0f1 = 0.0, f0 = 0.0, f1 = 0.0;
860 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
861 int t0 = GRAVI_BASE_TEL[base][0] + row_ft * ntel;
862 int t1 = GRAVI_BASE_TEL[base][1] + row_ft * ntel;
863 sf0f1 += sqrt (CPL_MAX(flux[t0] * flux[t1], 0));
869 if (f0==0 || f1==0)
continue;
872 pFactor[nsc] = sf0f1 * sf0f1 / (f0 * f1);
876 gravi_msg_function_exit(1);
877 return CPL_ERROR_NONE;
895 gravi_msg_function_start(1);
896 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
897 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
899 cpl_size nbase = 6, ntel = 4;
900 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) / ntel;
901 cpl_size nwave_sc = cpl_table_get_column_depth (flux_SC,
"FLUX");
904 cpl_array ** flux_sc = cpl_table_get_data_array (flux_SC,
"FLUX");
906 CPLCHECK_MSG (
"Cannot get data");
909 gravi_table_new_column_array (vis_SC,
"F1F2",
"e^2", CPL_TYPE_DOUBLE, nwave_sc);
910 cpl_array ** f1f2_sc = cpl_table_get_data_array (vis_SC,
"F1F2");
912 for (cpl_size base = 0; base < nbase; base++) {
913 for (cpl_size row = 0; row < nrow_sc; row ++) {
917 f1f2_sc[row*nbase+base] = cpl_array_cast (flux_sc[row*ntel+GRAVI_BASE_TEL[base][0]], CPL_TYPE_DOUBLE);
918 cpl_array_multiply (f1f2_sc[row*nbase+base], flux_sc[row*ntel+GRAVI_BASE_TEL[base][1]]);
922 gravi_msg_function_exit(1);
923 return CPL_ERROR_NONE;
943 gravi_msg_function_start(1);
944 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
945 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
947 cpl_size ntel = 4, nbase = 6;
948 cpl_size nrow_ft = cpl_table_get_nrow (flux_FT) / ntel;
949 cpl_size nwave_ft = cpl_table_get_column_depth (flux_FT,
"FLUX");
952 double * total_flux_ft = cpl_table_get_data_double (flux_FT,
"TOTALFLUX");
953 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
955 CPLCHECK_MSG (
"Cannot get data");
958 cpl_array ** mean_spectra = cpl_malloc (4 *
sizeof (cpl_array*));
960 for (cpl_size tel = 0; tel < ntel; tel ++) {
961 mean_spectra[tel] = cpl_array_duplicate (flux_ft[tel]);
962 for (cpl_size n = 1; n < nrow_ft; n ++) cpl_array_add (mean_spectra[tel], flux_ft[n*ntel+tel]);
963 cpl_array_divide_scalar (mean_spectra[tel], cpl_array_get_mean (mean_spectra[tel]) * nwave_ft);
964 CPLCHECK_MSG (
"Cannot compute mean spectra");
969 gravi_table_new_column_array (vis_FT,
"F1F2",
"e^2", CPL_TYPE_DOUBLE, nwave_ft);
970 cpl_array ** f1f2_ft = cpl_table_get_data_array (vis_FT,
"F1F2");
972 CPLCHECK_MSG (
"Cannot create columns");
975 for (cpl_size base = 0; base < nbase; base++) {
976 int t0 = GRAVI_BASE_TEL[base][0];
977 int t1 = GRAVI_BASE_TEL[base][1];
978 for (cpl_size n = 0; n < nrow_ft; n ++) {
979 f1f2_ft[n*nbase+base] = cpl_array_duplicate (mean_spectra[t0]);
980 cpl_array_multiply (f1f2_ft[n*nbase+base], mean_spectra[t1]);
981 cpl_array_multiply_scalar (f1f2_ft[n*nbase+base], total_flux_ft[n*ntel+t0] * total_flux_ft[n*ntel+t1]);
987 FREELOOP (cpl_array_delete, mean_spectra, 4);
989 gravi_msg_function_exit(1);
990 return CPL_ERROR_NONE;
1007 gravi_msg_function_start(1);
1008 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1012 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1013 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
1014 cpl_array ** visData_ft = cpl_table_get_data_array (vis_FT,
"VISDATA");
1016 CPLCHECK_MSG (
"Cannot get data");
1019 gravi_table_new_column_array (vis_FT,
"SELF_REF",
"rad", CPL_TYPE_DOUBLE, nwave_ft);
1020 cpl_array ** phaseref_ft = cpl_table_get_data_array (vis_FT,
"SELF_REF");
1022 CPLCHECK_MSG (
"Cannot create columns");
1026 for (cpl_size base = 0; base < nbase; base ++) {
1029 for (cpl_size n = 3; n < nrow_ft - 3; n ++) {
1030 phaseref_ft[n*nbase+base] = cpl_array_cast (visData_ft[(n-1)*nbase+base], CPL_TYPE_DOUBLE_COMPLEX);
1031 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+1)*nbase+base]);
1032 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], 2.0);
1033 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n-2)*nbase+base]);
1034 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+2)*nbase+base]);
1035 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], 2.0);
1036 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n-3)*nbase+base]);
1037 cpl_array_add (phaseref_ft[n*nbase+base], visData_ft[(n+3)*nbase+base]);
1039 cpl_array_arg (phaseref_ft[n*nbase+base]);
1040 cpl_array_multiply_scalar (phaseref_ft[n*nbase+base], -1.0);
1042 CPLCHECK_MSG(
"Cannot compute the PHASE_REF for the FT");
1044 phaseref_ft[0*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1045 phaseref_ft[1*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1046 phaseref_ft[2*nbase+base] = cpl_array_duplicate (phaseref_ft[3*nbase+base]);
1047 phaseref_ft[(nrow_ft-3)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1048 phaseref_ft[(nrow_ft-2)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1049 phaseref_ft[(nrow_ft-1)*nbase+base] = cpl_array_duplicate (phaseref_ft[(nrow_ft-4)*nbase+base]);
1053 cpl_msg_warning (cpl_func,
"Change of PHASE_REF_FT common channel new... to be decided !!");
1055 gravi_msg_function_exit(1);
1056 return CPL_ERROR_NONE;
1073 cpl_table * vis_ACQ)
1075 gravi_msg_function_start(1);
1076 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1077 cpl_ensure_code (vis_ACQ, CPL_ERROR_NULL_INPUT);
1079 cpl_size nbase = 6, ntel = 4;
1080 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1083 int * first = cpl_table_get_data_int (vis_SC,
"FIRST_ACQ");
1084 int * last = cpl_table_get_data_int (vis_SC,
"LAST_ACQ");
1085 CPLCHECK_MSG(
"Cannot get data");
1088 int * pup_n = cpl_table_get_data_int (vis_ACQ,
"PUPIL_NSPOT");
1089 double * pup_x = cpl_table_get_data_double (vis_ACQ,
"PUPIL_X");
1090 double * pup_y = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Y");
1091 double * pup_z = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Z");
1092 CPLCHECK_MSG(
"Cannot get direct pointer to data");
1095 gravi_table_new_column (vis_SC,
"PUPIL_X", NULL, CPL_TYPE_DOUBLE);
1096 double * pup_x_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_X");
1097 gravi_table_new_column (vis_SC,
"PUPIL_Y", NULL, CPL_TYPE_DOUBLE);
1098 double * pup_y_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_Y");
1099 gravi_table_new_column (vis_SC,
"PUPIL_Z", NULL, CPL_TYPE_DOUBLE);
1100 double * pup_z_sc = cpl_table_get_data_double (vis_SC,
"PUPIL_Z");
1103 for (cpl_size base = 0; base < nbase; base++) {
1104 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1105 cpl_size nsc = row_sc * nbase + base;
1108 cpl_size nframe = 0;
1109 for (cpl_size row = first[nsc] ; row < last[nsc]; row++) {
1110 cpl_size row0 = row * ntel + GRAVI_BASE_TEL[base][0];
1111 cpl_size row1 = row * ntel + GRAVI_BASE_TEL[base][1];
1113 if (pup_n[row0] != 0 && pup_n[row1] !=0 ) {
1114 pup_x_sc[nsc] += pup_x[row0] - pup_x[row1];
1115 pup_y_sc[nsc] += pup_y[row0] - pup_y[row1];
1116 pup_z_sc[nsc] += pup_z[row0] - pup_z[row1];
1120 CPLCHECK_MSG (
"Fail to integrate the ACQ frames");
1125 pup_x_sc[nsc] /= (double)nframe;
1126 pup_y_sc[nsc] /= (double)nframe;
1127 pup_z_sc[nsc] /= (double)nframe;
1133 gravi_msg_function_exit(1);
1134 return CPL_ERROR_NONE;
1151 cpl_table * vis_ACQ)
1153 gravi_msg_function_start(1);
1154 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1155 cpl_ensure_code (vis_ACQ, CPL_ERROR_NULL_INPUT);
1158 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) / ntel;
1161 int * first = cpl_table_get_data_int (flux_SC,
"FIRST_ACQ");
1162 int * last = cpl_table_get_data_int (flux_SC,
"LAST_ACQ");
1163 CPLCHECK_MSG (
"Cannot get data");
1166 int * pup_n = cpl_table_get_data_int (vis_ACQ,
"PUPIL_NSPOT");
1167 double * pup_x = cpl_table_get_data_double (vis_ACQ,
"PUPIL_X");
1168 double * pup_y = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Y");
1169 double * pup_z = cpl_table_get_data_double (vis_ACQ,
"PUPIL_Z");
1170 CPLCHECK_MSG (
"Cannot get direct pointer to data");
1173 gravi_table_new_column (flux_SC,
"PUPIL_X", NULL, CPL_TYPE_DOUBLE);
1174 double * pup_x_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_X");
1175 gravi_table_new_column (flux_SC,
"PUPIL_Y", NULL, CPL_TYPE_DOUBLE);
1176 double * pup_y_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_Y");
1177 gravi_table_new_column (flux_SC,
"PUPIL_Z", NULL, CPL_TYPE_DOUBLE);
1178 double * pup_z_sc = cpl_table_get_data_double (flux_SC,
"PUPIL_Z");
1181 for (cpl_size tel = 0; tel < ntel; tel++) {
1182 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1183 cpl_size nsc = row_sc * ntel + tel;
1186 cpl_size nframe = 0;
1187 for (cpl_size row = first[nsc] ; row < last[nsc]; row++) {
1188 cpl_size row0 = row * ntel + tel;
1190 if (pup_n[row0] != 0) {
1191 pup_x_sc[nsc] += pup_x[row0];
1192 pup_y_sc[nsc] += pup_y[row0];
1193 pup_z_sc[nsc] += pup_z[row0];
1197 CPLCHECK_MSG (
"Fail to integrate the ACQ frames");
1202 pup_x_sc[nsc] /= (double)nframe;
1203 pup_y_sc[nsc] /= (double)nframe;
1204 pup_z_sc[nsc] /= (double)nframe;
1210 gravi_msg_function_exit(1);
1211 return CPL_ERROR_NONE;
1229 gravi_msg_function_start(1);
1230 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1231 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1233 cpl_size nbase = 6, ndiode = 4, ntel = 4;
1234 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1237 int * first_met = cpl_table_get_data_int (vis_SC,
"FIRST_MET");
1238 int * last_met = cpl_table_get_data_int (vis_SC,
"LAST_MET");
1240 CPLCHECK_MSG(
"Cannot get data");
1243 double * phase_met_fc = cpl_table_get_data_double (vis_MET,
"PHASE_FC");
1244 double * opd_met_fc = cpl_table_get_data_double (vis_MET,
"OPD_FC");
1245 cpl_array ** phase_met_tel = cpl_table_get_data_array (vis_MET,
"PHASE_TEL");
1246 cpl_array ** opd_met_tel = cpl_table_get_data_array (vis_MET,
"OPD_TEL");
1247 cpl_array ** phasor_met_telfc = cpl_table_get_data_array (vis_MET,
"PHASOR_TELFC");
1249 double * opd_met_fc_corr = cpl_table_get_data_double (vis_MET,
"OPD_FC_CORR");
1250 double * opd_met_telfc_mcorr = cpl_table_get_data_double (vis_MET,
"OPD_TELFC_MCORR");
1251 cpl_array ** opd_met_telfc_corr = cpl_table_get_data_array (vis_MET,
"OPD_TELFC_CORR");
1253 CPLCHECK_MSG(
"Cannot get direct pointer to data");
1256 gravi_table_new_column (vis_SC,
"PHASE_MET_FC",
"rad", CPL_TYPE_DOUBLE);
1257 double * phase_metdit_fc = cpl_table_get_data_double (vis_SC,
"PHASE_MET_FC");
1259 gravi_table_new_column (vis_SC,
"PHASE_MET_TEL",
"rad", CPL_TYPE_DOUBLE);
1260 double * phase_metdit_tel = cpl_table_get_data_double (vis_SC,
"PHASE_MET_TEL");
1262 gravi_table_new_column (vis_SC,
"OPD_MET_FC",
"m", CPL_TYPE_DOUBLE);
1263 double * opd_metdit_fc = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC");
1265 gravi_table_new_column_array (vis_SC,
"OPD_MET_TEL",
"m", CPL_TYPE_DOUBLE, ndiode);
1266 cpl_array ** opd_metdit_tel = cpl_table_get_data_array (vis_SC,
"OPD_MET_TEL");
1268 gravi_table_new_column_array (vis_SC,
"PHASOR_MET_TELFC",
"V^4", CPL_TYPE_DOUBLE_COMPLEX, ndiode);
1269 cpl_array ** phasor_metdit_telfc = cpl_table_get_data_array (vis_SC,
"PHASOR_MET_TELFC");
1271 gravi_table_new_column (vis_SC,
"OPD_MET_FC_CORR",
"m", CPL_TYPE_DOUBLE);
1272 double * opd_metdit_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
1274 gravi_table_new_column (vis_SC,
"OPD_MET_TELFC_MCORR",
"m", CPL_TYPE_DOUBLE);
1275 double * opd_metdit_telfc_mcorr = cpl_table_get_data_double (vis_SC,
"OPD_MET_TELFC_MCORR");
1277 gravi_table_new_column_array (vis_SC,
"OPD_MET_TELFC_CORR",
"m", CPL_TYPE_DOUBLE, ndiode);
1278 cpl_array ** opd_metdit_telfc_corr = cpl_table_get_data_array (vis_SC,
"OPD_MET_TELFC_CORR");
1280 CPLCHECK_MSG(
"Cannot create columns");
1283 for (cpl_size base = 0; base < nbase; base++) {
1284 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1285 cpl_size nsc = row_sc * nbase + base;
1287 opd_metdit_tel[nsc] = gravi_array_init_double (ndiode, 0.0);
1288 phasor_metdit_telfc[nsc] = gravi_array_init_double_complex (ndiode, 0.0+I*0.0);
1289 opd_metdit_telfc_corr[nsc] = gravi_array_init_double (ndiode, 0.0);
1292 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1293 cpl_size nmet0 = row_met * ntel + GRAVI_BASE_TEL[base][0];
1294 cpl_size nmet1 = row_met * ntel + GRAVI_BASE_TEL[base][1];
1297 opd_metdit_fc_corr[nsc] += opd_met_fc_corr[nmet0] - opd_met_fc_corr[nmet1];
1298 opd_metdit_telfc_mcorr[nsc] += opd_met_telfc_mcorr[nmet0] - opd_met_telfc_mcorr[nmet1];
1301 cpl_array_add (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet0]);
1302 cpl_array_subtract (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet1]);
1305 cpl_array_add (opd_metdit_tel[nsc], opd_met_tel[nmet0]);
1306 cpl_array_subtract (opd_metdit_tel[nsc], opd_met_tel[nmet1]);
1309 phase_metdit_tel[nsc] += cpl_array_get_mean (phase_met_tel[nmet0]) -
1310 cpl_array_get_mean (phase_met_tel[nmet1]);
1313 phase_metdit_fc[nsc] += phase_met_fc[nmet0] - phase_met_fc[nmet1];
1316 opd_metdit_fc[nsc] += opd_met_fc[nmet0] - opd_met_fc[nmet1];
1320 phasor_met_telfc[nmet0],
1321 phasor_met_telfc[nmet1]);
1323 CPLCHECK_MSG (
"Fail to integrate the metrology");
1327 cpl_size nframe = last_met[nsc] - first_met[nsc];
1329 opd_metdit_fc_corr[nsc] /= nframe;
1330 opd_metdit_telfc_mcorr[nsc] /= nframe;
1331 cpl_array_divide_scalar (opd_metdit_telfc_corr[nsc], (
double)nframe);
1333 cpl_array_divide_scalar (opd_metdit_tel[nsc], (
double)nframe);
1334 phase_metdit_tel[nsc] /= nframe;
1335 phase_metdit_fc[nsc] /= nframe;
1336 opd_metdit_fc[nsc] /= nframe;
1337 cpl_array_divide_scalar (phasor_metdit_telfc[nsc], (
double)nframe);
1339 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1347 if (cpl_table_has_column (vis_MET,
"FIELD_FIBER_DX")) {
1349 double * fdx_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DX");
1350 double * fdy_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DY");
1352 gravi_table_new_column (vis_SC,
"FIELD_FIBER_DX",
"pix", CPL_TYPE_DOUBLE);
1353 double * fdx_metdit = cpl_table_get_data_double (vis_SC,
"FIELD_FIBER_DX");
1355 gravi_table_new_column (vis_SC,
"FIELD_FIBER_DY",
"pix", CPL_TYPE_DOUBLE);
1356 double * fdy_metdit = cpl_table_get_data_double (vis_SC,
"FIELD_FIBER_DY");
1359 for (cpl_size base = 0; base < nbase; base++) {
1360 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1361 cpl_size nsc = row_sc * nbase + base;
1364 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1365 cpl_size nmet0 = row_met * ntel + GRAVI_BASE_TEL[base][0];
1366 cpl_size nmet1 = row_met * ntel + GRAVI_BASE_TEL[base][1];
1369 fdx_metdit[nsc] += fdx_met[nmet0] - fdx_met[nmet1];
1370 fdy_metdit[nsc] += fdy_met[nmet0] - fdy_met[nmet1];
1372 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1375 cpl_size nframe = last_met[nsc] - first_met[nsc];
1377 fdx_metdit[nsc] /= nframe;
1378 fdy_metdit[nsc] /= nframe;
1386 gravi_msg_function_exit(1);
1387 return CPL_ERROR_NONE;
1414 gravi_msg_function_start(1);
1415 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1416 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1418 cpl_size nbase = 6, ntel = 4;
1419 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1420 cpl_size nrow_met = cpl_table_get_nrow (vis_MET) / ntel;
1422 CPLCHECK_MSG(
"Cannot get data");
1425 double periode_ft = cpl_table_get (vis_FT,
"TIME", nbase, NULL) -
1426 cpl_table_get (vis_FT,
"TIME", 0, NULL);
1427 cpl_msg_info (cpl_func,
"PERIOD FT = %g [us]", periode_ft);
1430 double * phase_met_fc = cpl_table_get_data_double (vis_MET,
"PHASE_FC");
1432 CPLCHECK_MSG(
"Cannot get direct pointer to data");
1435 gravi_table_new_column (vis_FT,
"PHASE_MET_FC",
"rad", CPL_TYPE_DOUBLE);
1436 double * phase_metdit_fc = cpl_table_get_data_double (vis_FT,
"PHASE_MET_FC");
1438 CPLCHECK_MSG(
"Cannot create columns");
1441 for (cpl_size base = 0; base < nbase; base++) {
1442 int tel0 = GRAVI_BASE_TEL[base][0];
1443 int tel1 = GRAVI_BASE_TEL[base][1];
1445 for (cpl_size last_met = 0, first_met = 0, row_ft = 0; row_ft < nrow_ft; row_ft ++) {
1446 cpl_size nft = row_ft * nbase + base;
1447 double time_ft = cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL);
1456 first_met = CPL_MAX (CPL_MIN (last_met - 5, nrow_met - 1), 0);
1457 while ((cpl_table_get (vis_MET,
"TIME", first_met*ntel, NULL) <= (time_ft + 0.0))) {
1459 if (first_met == nrow_met)
break;
1461 CPLCHECK_MSG (
"Cannot get first");
1464 last_met = CPL_MAX (CPL_MIN (first_met - 1, nrow_met - 1), 0);
1465 while ((cpl_table_get (vis_MET,
"TIME", last_met*ntel, NULL) < (time_ft + periode_ft))) {
1467 if (last_met == nrow_met)
break;
1469 CPLCHECK_MSG (
"Cannot get last");
1472 if (row_ft < 5 && last_met == 0) last_met = 1;
1476 if (row_ft > nrow_ft-5 && first_met == nrow_met) first_met = nrow_met - 1;
1477 if (row_ft > nrow_ft-5 && last_met == nrow_met) last_met = nrow_met;
1480 if ( last_met == 0 || last_met > nrow_met ) {
1481 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
1482 "Not enough MET data to synchronise " 1483 "with FT DIT %lli over %lli", row_ft+1, nrow_ft);
1486 if ( last_met - first_met > 0 ) {
1490 for (cpl_size rin_met = first_met; rin_met < last_met; rin_met ++)
1491 phase_metdit_fc[nft] += phase_met_fc[rin_met*ntel+tel0] - phase_met_fc[rin_met*ntel+tel1];
1493 phase_metdit_fc[nft] /= (last_met - first_met);
1497 cpl_size rowa_met = first_met-1;
1498 cpl_size rowb_met = last_met;
1500 double phia_met = phase_met_fc[rowa_met*ntel+tel0] - phase_met_fc[rowa_met*ntel+tel1];
1501 double phib_met = phase_met_fc[rowb_met*ntel+tel0] - phase_met_fc[rowb_met*ntel+tel1];
1502 double timea_met = cpl_table_get (vis_MET,
"TIME",rowa_met*ntel, NULL);
1503 double timeb_met = cpl_table_get (vis_MET,
"TIME",rowb_met*ntel, NULL);
1505 phase_metdit_fc[nft] = phia_met +
1506 (phib_met - phia_met) * (time_ft - timea_met) / (timeb_met - timea_met);
1507 CPLCHECK_MSG (
"Cannot interpolate");
1510 CPLCHECK_MSG (
"Fail to compute metrology per FT base from metrology per tel");
1515 gravi_msg_function_exit(1);
1516 return CPL_ERROR_NONE;
1542 gravi_msg_function_start(1);
1543 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1544 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1547 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
1548 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1550 CPLCHECK_MSG(
"Cannot get data");
1553 double * phase_sc = cpl_table_get_data_double (vis_SC,
"OPD");
1555 CPLCHECK_MSG(
"Cannot get direct pointer to data");
1558 gravi_table_new_column (vis_FT,
"OPD_SC",
"rad", CPL_TYPE_DOUBLE);
1559 double * phase_scdit = cpl_table_get_data_double (vis_FT,
"OPD_SC");
1561 CPLCHECK_MSG(
"Cannot create columns");
1564 for (cpl_size base = 0; base < nbase; base++) {
1565 cpl_size row_ft = 0;
1567 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1568 cpl_size nsc = row_sc * nbase + base;
1569 double time_sc = cpl_table_get (vis_SC,
"TIME", nsc, NULL);
1572 while ((cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL) < (time_sc - dit_sc/2))) {
1574 if (row_ft >= nrow_ft)
break;
1578 while ((cpl_table_get (vis_FT,
"TIME", row_ft*nbase+base, NULL) < (time_sc + dit_sc/2))) {
1580 phase_scdit[row_ft*nbase+base] = phase_sc[nsc];
1583 if (row_ft >= nrow_ft)
break;
1586 if (row_ft >= nrow_ft) {
1587 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
1588 "Not enough FT data to synchronise " 1589 "with SC DIT %lli over %lli", row_sc+1, nrow_sc);
1596 gravi_msg_function_exit(1);
1597 return CPL_ERROR_NONE;
1616 gravi_msg_function_start(1);
1617 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1618 cpl_ensure_code (vis_MET, CPL_ERROR_NULL_INPUT);
1620 cpl_size ntel = 4, ndiode = 4;
1621 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) / ntel;
1624 int * first_met = cpl_table_get_data_int (flux_SC,
"FIRST_MET");
1625 int * last_met = cpl_table_get_data_int (flux_SC,
"LAST_MET");
1627 CPLCHECK_MSG(
"Cannot get data");
1630 double * opd_met_fc = cpl_table_get_data_double (vis_MET,
"OPD_FC");
1631 cpl_array ** opd_met_tel = cpl_table_get_data_array (vis_MET,
"OPD_TEL");
1632 cpl_array ** phasor_met_telfc = cpl_table_get_data_array (vis_MET,
"PHASOR_TELFC");
1634 double * opd_met_fc_corr = cpl_table_get_data_double (vis_MET,
"OPD_FC_CORR");
1635 double * opd_met_telfc_mcorr = cpl_table_get_data_double (vis_MET,
"OPD_TELFC_MCORR");
1636 cpl_array ** opd_met_telfc_corr = cpl_table_get_data_array (vis_MET,
"OPD_TELFC_CORR");
1638 CPLCHECK_MSG(
"Cannot get direct pointer to data");
1641 gravi_table_new_column (flux_SC,
"OPD_MET_FC",
"m", CPL_TYPE_DOUBLE);
1642 double * opd_metdit_fc = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC");
1644 gravi_table_new_column_array (flux_SC,
"OPD_MET_TEL",
"m", CPL_TYPE_DOUBLE, ndiode);
1645 cpl_array ** opd_metdit_tel = cpl_table_get_data_array (flux_SC,
"OPD_MET_TEL");
1647 gravi_table_new_column_array (flux_SC,
"PHASOR_MET_TELFC",
"V^4", CPL_TYPE_DOUBLE_COMPLEX, ndiode);
1648 cpl_array ** phasor_metdit_telfc = cpl_table_get_data_array (flux_SC,
"PHASOR_MET_TELFC");
1650 gravi_table_new_column (flux_SC,
"OPD_MET_FC_CORR",
"m", CPL_TYPE_DOUBLE);
1651 double * opd_metdit_fc_corr = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC_CORR");
1653 gravi_table_new_column (flux_SC,
"OPD_MET_TELFC_MCORR",
"m", CPL_TYPE_DOUBLE);
1654 double * opd_metdit_telfc_mcorr = cpl_table_get_data_double (flux_SC,
"OPD_MET_TELFC_MCORR");
1656 gravi_table_new_column_array (flux_SC,
"OPD_MET_TELFC_CORR",
"m", CPL_TYPE_DOUBLE, ndiode);
1657 cpl_array ** opd_metdit_telfc_corr = cpl_table_get_data_array (flux_SC,
"OPD_MET_TELFC_CORR");
1659 CPLCHECK_MSG(
"Cannot create columns");
1662 for (cpl_size tel = 0; tel < ntel; tel++) {
1663 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1664 cpl_size nsc = row_sc * ntel + tel;
1666 opd_metdit_tel[nsc] = gravi_array_init_double (ndiode, 0.0);
1667 phasor_metdit_telfc[nsc] = gravi_array_init_double_complex (ndiode, 0.0+I*0.0);
1668 opd_metdit_telfc_corr[nsc] = gravi_array_init_double (ndiode, 0.0);
1671 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1672 cpl_size nmet = row_met * ntel + tel;
1675 cpl_array_add (opd_metdit_tel[nsc], opd_met_tel[nmet]);
1678 opd_metdit_fc[nsc] += opd_met_fc[nmet];
1681 cpl_array_add (phasor_metdit_telfc[nsc], phasor_met_telfc[nmet]);
1684 opd_metdit_fc_corr[nsc] += opd_met_fc_corr[nmet];
1685 opd_metdit_telfc_mcorr[nsc] += opd_met_telfc_mcorr[nmet];
1688 cpl_array_add (opd_metdit_telfc_corr[nsc], opd_met_telfc_corr[nmet]);
1690 CPLCHECK_MSG (
"Fail to integrate the metrology");
1694 cpl_size nframe = last_met[nsc] - first_met[nsc];
1696 opd_metdit_fc_corr[nsc] /= nframe;
1697 opd_metdit_telfc_mcorr[nsc] /= nframe;
1698 cpl_array_divide_scalar (opd_metdit_telfc_corr[nsc], (
double)nframe);
1700 cpl_array_divide_scalar (opd_metdit_tel[nsc], (
double)nframe);
1701 opd_metdit_fc[nsc] /= nframe;
1702 cpl_array_divide_scalar (phasor_metdit_telfc[nsc], (
double)nframe);
1704 CPLCHECK_MSG (
"Fail to integrate the metrology");
1712 if (cpl_table_has_column (vis_MET,
"FIELD_FIBER_DX")) {
1714 double * fdx_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DX");
1715 double * fdy_met = cpl_table_get_data_double (vis_MET,
"FIELD_FIBER_DY");
1717 gravi_table_new_column (flux_SC,
"FIELD_FIBER_DX",
"pix", CPL_TYPE_DOUBLE);
1718 double * fdx_metdit = cpl_table_get_data_double (flux_SC,
"FIELD_FIBER_DX");
1720 gravi_table_new_column (flux_SC,
"FIELD_FIBER_DY",
"pix", CPL_TYPE_DOUBLE);
1721 double * fdy_metdit = cpl_table_get_data_double (flux_SC,
"FIELD_FIBER_DY");
1724 for (cpl_size tel = 0; tel < ntel; tel++) {
1725 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1726 cpl_size nsc = row_sc * ntel + tel;
1729 for (cpl_size row_met = first_met[nsc] ; row_met < last_met[nsc]; row_met++) {
1730 cpl_size nmet0 = row_met * ntel + tel;
1733 fdx_metdit[nsc] += fdx_met[nmet0];
1734 fdy_metdit[nsc] += fdy_met[nmet0];
1736 CPLCHECK_MSG (
"Fail to compute metrology per base from metrology per tel");
1739 cpl_size nframe = last_met[nsc] - first_met[nsc];
1741 fdx_metdit[nsc] /= nframe;
1742 fdy_metdit[nsc] /= nframe;
1749 gravi_msg_function_exit(1);
1750 return CPL_ERROR_NONE;
1769 gravi_msg_function_start(1);
1770 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1771 cpl_ensure_code (fddl_table, CPL_ERROR_NULL_INPUT);
1774 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) / ntel;
1777 int * first_fddl = cpl_table_get_data_int (flux_SC,
"FIRST_FDDL");
1778 int * last_fddl = cpl_table_get_data_int (flux_SC,
"LAST_FDDL");
1780 cpl_array ** ftpos = cpl_table_get_data_array (fddl_table,
"FT_POS");
1781 cpl_array ** scpos = cpl_table_get_data_array (fddl_table,
"SC_POS");
1782 cpl_array ** oplair = cpl_table_get_data_array (fddl_table,
"OPL_AIR");
1784 CPLCHECK_MSG(
"Cannot get data");
1787 gravi_table_new_column (flux_SC,
"FT_POS",
"V", CPL_TYPE_DOUBLE);
1788 double * ftpos_flux = cpl_table_get_data_double (flux_SC,
"FT_POS");
1790 gravi_table_new_column (flux_SC,
"SC_POS",
"V", CPL_TYPE_DOUBLE);
1791 double * scpos_flux = cpl_table_get_data_double (flux_SC,
"SC_POS");
1793 gravi_table_new_column (flux_SC,
"OPL_AIR",
"m", CPL_TYPE_DOUBLE);
1794 double * oplair_flux = cpl_table_get_data_double (flux_SC,
"OPL_AIR");
1796 CPLCHECK_MSG(
"Cannot create columns");
1799 for (cpl_size tel = 0; tel < ntel; tel++) {
1800 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1801 cpl_size nsc = row_sc * ntel + tel;
1804 for (cpl_size row_fddl = first_fddl[nsc] ; row_fddl < last_fddl[nsc]; row_fddl++) {
1805 ftpos_flux[nsc] += cpl_array_get (ftpos[row_fddl], tel, NULL);
1806 scpos_flux[nsc] += cpl_array_get (scpos[row_fddl], tel, NULL);
1807 oplair_flux[nsc] += cpl_array_get (oplair[row_fddl], tel, NULL);
1811 cpl_size nframe = last_fddl[nsc] - first_fddl[nsc];
1813 ftpos_flux[nsc] /= nframe;
1814 scpos_flux[nsc] /= nframe;
1815 oplair_flux[nsc] /= nframe;
1821 gravi_msg_function_exit(1);
1822 return CPL_ERROR_NONE;
1840 gravi_msg_function_start(1);
1841 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
1842 cpl_ensure_code (flux_FT, CPL_ERROR_NULL_INPUT);
1845 cpl_size nrow_sc = cpl_table_get_nrow (flux_SC) / ntel;
1846 cpl_size nwave_sc = cpl_table_get_column_depth (flux_SC,
"FLUX");
1847 cpl_size nwave_ft = cpl_table_get_column_depth (flux_FT,
"FLUX");
1850 int * first_ft = cpl_table_get_data_int (flux_SC,
"FIRST_FT");
1851 int * last_ft = cpl_table_get_data_int (flux_SC,
"LAST_FT");
1853 cpl_array ** flux_sc = cpl_table_get_data_array (flux_SC,
"FLUX");
1854 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
1856 CPLCHECK_MSG(
"Cannot get data");
1859 gravi_table_new_column (flux_SC,
"TOTALFLUX_SC",
"e", CPL_TYPE_DOUBLE);
1860 double * total_flux_scdit = cpl_table_get_data_double (flux_SC,
"TOTALFLUX_SC");
1862 gravi_table_new_column (flux_SC,
"TOTALFLUX_FT",
"e", CPL_TYPE_DOUBLE);
1863 double * total_flux_ftdit = cpl_table_get_data_double (flux_SC,
"TOTALFLUX_FT");
1865 CPLCHECK_MSG(
"Cannot create columns");
1867 for (cpl_size tel = 0; tel < ntel; tel++) {
1868 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
1869 cpl_size nsc = row_sc * ntel + tel;
1872 total_flux_scdit[nsc] = cpl_array_get_mean (flux_sc[nsc]) * nwave_sc;
1875 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
1876 total_flux_ftdit[nsc] += cpl_array_get_mean (flux_ft[row_ft * ntel + tel]) * nwave_ft;
1879 CPLCHECK_MSG(
"Issue in the loop to average FT and SC flux per frame");
1884 gravi_msg_function_exit(1);
1885 return CPL_ERROR_NONE;
1907 cpl_table * wave_table_sc,
1909 cpl_table * wave_table_ft)
1911 gravi_msg_function_start(1);
1912 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
1913 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
1914 cpl_ensure_code (wave_table_sc, CPL_ERROR_NULL_INPUT);
1915 cpl_ensure_code (wave_table_ft, CPL_ERROR_NULL_INPUT);
1919 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
1920 cpl_size nwave_sc = cpl_table_get_column_depth (vis_SC,
"VISDATA");
1921 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
1924 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
1925 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
1927 cpl_array ** visErr_ft = cpl_table_get_data_array (vis_FT,
"VISERR");
1928 cpl_array ** visData_ft = cpl_table_get_data_array (vis_FT,
"VISDATA");
1930 CPLCHECK_MSG (
"Cannot get data");
1933 double meanwave_ft = cpl_table_get_column_mean (wave_table_ft,
"EFF_WAVE");
1934 double * wavenumber_sc = cpl_malloc (nwave_sc *
sizeof(
double));
1935 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
1936 wavenumber_sc[wave] = 1. / cpl_table_get (wave_table_sc,
"EFF_WAVE", wave, NULL);
1939 CPLCHECK_MSG (
"Cannot get data");
1942 gravi_table_new_column_array (vis_SC,
"VISDATA_FT",
"e", CPL_TYPE_DOUBLE_COMPLEX, nwave_ft);
1943 cpl_array ** visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA_FT");
1945 gravi_table_new_column_array (vis_SC,
"VISVAR_FT",
"e^2", CPL_TYPE_DOUBLE, nwave_ft);
1946 cpl_array ** visVar_ftdit = cpl_table_get_data_array (vis_SC,
"VISVAR_FT");
1948 gravi_table_new_column_array (vis_SC,
"VISPOWER_FT",
"e^2", CPL_TYPE_DOUBLE, nwave_ft);
1949 cpl_array ** visPower_ftdit = cpl_table_get_data_array (vis_SC,
"VISPOWER_FT");
1951 gravi_table_new_column_array (vis_SC,
"V_FACTOR", NULL, CPL_TYPE_DOUBLE, nwave_sc);
1952 cpl_array ** vFactor = cpl_table_get_data_array (vis_SC,
"V_FACTOR");
1954 gravi_table_new_column_array (vis_SC,
"V_FACTOR_FT", NULL, CPL_TYPE_DOUBLE, nwave_ft);
1955 cpl_array ** vFactor_ftdit = cpl_table_get_data_array (vis_SC,
"V_FACTOR_FT");
1957 gravi_table_new_column (vis_SC,
"V_FACTOR_WL", NULL, CPL_TYPE_DOUBLE);
1958 double * vFactor_wl = cpl_table_get_data_double (vis_SC,
"V_FACTOR_WL");
1960 CPLCHECK_MSG (
"Cannot create columns");
1963 for (cpl_size base = 0; base < nbase; base++) {
1964 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
1965 int nsc = row_sc * nbase + base;
1968 visData_ftdit[nsc] = gravi_array_init_double_complex (nwave_ft, 0.0 + I*0.0);
1969 visVar_ftdit[nsc] = gravi_array_init_double (nwave_ft, 0.0);
1970 visPower_ftdit[nsc] = gravi_array_init_double (nwave_ft, 0.0);
1972 double vFactor_incoh = 0.0;
1973 double complex vFactor_coh = 0.0 + I * 0.0;
1974 double vFactor_var = 0.0;
1977 cpl_size nframe = last_ft[nsc] - first_ft[nsc];
1978 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
1979 cpl_array * tmp_data;
1982 tmp_data = visErr_ft[row_ft * nbase + base];
1983 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
1984 cpl_array_set (visVar_ftdit[nsc], wave,
1985 cpl_array_get (visVar_ftdit[nsc], wave, NULL) +
1986 pow (cabs( cpl_array_get_complex (tmp_data, wave, NULL) ), 2));
1990 tmp_data = visData_ft[row_ft * nbase + base];
1991 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
1992 cpl_array_set_complex (visData_ftdit[nsc], wave,
1993 cpl_array_get_complex (visData_ftdit[nsc], wave, NULL) +
1994 cpl_array_get_complex (tmp_data, wave, NULL));
1998 for (cpl_size wave = 0; wave < nwave_ft; wave ++){
1999 cpl_array_set (visPower_ftdit[nsc], wave,
2000 cpl_array_get (visPower_ftdit[nsc], wave, NULL) +
2001 pow (cabs( cpl_array_get_complex (tmp_data, wave, NULL)), 2));
2006 vFactor_incoh += pow (cabs (cpl_array_get_mean_complex (visData_ft[row_ft * nbase + base])) * nwave_ft, 2);
2007 vFactor_coh += cpl_array_get_mean_complex (visData_ft[row_ft * 6 + base]) * nwave_ft;
2008 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2009 vFactor_var += pow (cabs (cpl_array_get_complex (visErr_ft[row_ft * nbase + base], wave, NULL)), 2);
2012 CPLCHECK_MSG(
"Issue in the loop to build average FT quantities");
2017 vFactor_ftdit[nsc] = cpl_array_new (nwave_ft, CPL_TYPE_DOUBLE);
2019 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2020 cpl_array_set (vFactor_ftdit[nsc], wave,
2021 (pow (cabs (cpl_array_get_complex (visData_ftdit[nsc], wave, NULL)),2) -
2022 cpl_array_get (visVar_ftdit[nsc], wave, NULL)) /
2023 (cpl_array_get (visPower_ftdit[nsc], wave, NULL) -
2024 cpl_array_get (visVar_ftdit[nsc], wave, NULL)) / (
double)nframe);
2029 vFactor_wl[nsc] = (cabs(vFactor_coh)*cabs(vFactor_coh) - vFactor_var) /
2030 (vFactor_incoh - vFactor_var) / (double)nframe;
2034 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2035 cpl_array_set (vFactor_ftdit[nsc], wave,0);
2037 vFactor_wl[nsc] = 0;
2042 vFactor[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2043 double a2 = log (CPL_MAX(CPL_MIN (vFactor_wl[nsc], 1.0),1e-10)) * meanwave_ft*meanwave_ft;
2044 for (cpl_size wave = 0; wave < nwave_sc; wave ++) {
2045 cpl_array_set (vFactor[nsc], wave, exp ( a2 * pow (wavenumber_sc[wave], 2)));
2052 FREE (cpl_free, wavenumber_sc);
2054 gravi_msg_function_exit(1);
2055 return CPL_ERROR_NONE;
2077 gravi_msg_function_start(1);
2078 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2079 cpl_ensure_code (vis_FT, CPL_ERROR_NULL_INPUT);
2083 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2086 int * reject_flag_ft = cpl_table_get_data_int (vis_FT,
"REJECTION_FLAG");
2087 int * first_ft = cpl_table_get_data_int (vis_SC,
"FIRST_FT");
2088 int * last_ft = cpl_table_get_data_int (vis_SC,
"LAST_FT");
2090 CPLCHECK_MSG (
"Cannot get data");
2093 gravi_table_new_column (vis_SC,
"FRINGEDET_RATIO", NULL, CPL_TYPE_DOUBLE);
2094 double * fringedet_ftdit = cpl_table_get_data_double (vis_SC,
"FRINGEDET_RATIO");
2096 CPLCHECK_MSG (
"Cannot create columns");
2099 for (cpl_size base = 0; base < nbase; base++) {
2100 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc++) {
2101 int nsc = row_sc * nbase + base;
2104 for (cpl_size row_ft = first_ft[nsc] ; row_ft < last_ft[nsc]; row_ft++) {
2105 fringedet_ftdit[nsc] += (reject_flag_ft[row_ft * nbase + base] == 0 ? 1 : 0);
2109 cpl_size nframe = last_ft[nsc] - first_ft[nsc];
2110 if (nframe != 0) fringedet_ftdit[nsc] /= (double)nframe;
2114 gravi_msg_function_exit(1);
2115 return CPL_ERROR_NONE;
2138 cpl_table * wavesc_table,
2139 cpl_table * waveft_table)
2141 gravi_msg_function_start(1);
2142 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2143 cpl_ensure_code (wavesc_table, CPL_ERROR_NULL_INPUT);
2145 cpl_array ** visData_ftdit;
2146 double * gdelay_ftdit;
2147 const char * output_name = NULL, * coeff_name = NULL;
2148 if (waveft_table != NULL) {
2149 cpl_msg_info (cpl_func,
"Compute reference phase of SC from the FT data");
2150 visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA_FT");
2151 gdelay_ftdit = cpl_table_get_data_double (vis_SC,
"GDELAY_FT");
2152 output_name =
"PHASE_REF";
2153 coeff_name =
"PHASE_REF_COEFF";
2155 cpl_msg_info (cpl_func,
"Compute reference phase of SC from the SC");
2156 visData_ftdit = cpl_table_get_data_array (vis_SC,
"VISDATA");
2157 gdelay_ftdit = cpl_table_get_data_double (vis_SC,
"GDELAY");
2158 waveft_table = wavesc_table;
2159 output_name =
"SELF_REF";
2160 coeff_name =
"SELF_REF_COEFF";
2165 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2166 cpl_size nwave_sc = cpl_table_get_nrow (wavesc_table);
2167 cpl_size nwave_ft = cpl_table_get_nrow (waveft_table);
2169 CPLCHECK_MSG (
"Cannot get data");
2172 cpl_size mindeg = 0, maxdeg = 2;
2173 cpl_polynomial * fit = cpl_polynomial_new (1);
2176 cpl_matrix * wave_ft = cpl_matrix_new (1,nwave_ft);
2177 cpl_vector * wave_sc = cpl_vector_new (nwave_sc);
2178 cpl_array * wavenumber_ft = cpl_array_new (nwave_ft, CPL_TYPE_DOUBLE);
2180 double lbd0 = cpl_table_get_column_mean (waveft_table,
"EFF_WAVE");
2181 double delta0 = cpl_table_get_column_max (waveft_table,
"EFF_WAVE") -
2182 cpl_table_get_column_min (waveft_table,
"EFF_WAVE");
2183 for (cpl_size wave = 0; wave < nwave_ft; wave ++) {
2184 double lbd = cpl_table_get (waveft_table,
"EFF_WAVE", wave, NULL);
2185 cpl_matrix_set (wave_ft, 0, wave, (lbd - lbd0) / delta0);
2186 cpl_array_set (wavenumber_ft, wave, 1./lbd);
2188 for (cpl_size wave = 0; wave < nwave_sc; wave ++) {
2189 cpl_vector_set (wave_sc, wave, cpl_table_get (wavesc_table,
"EFF_WAVE", wave, NULL));
2192 CPLCHECK_MSG (
"Cannot create wave arrays");
2195 gravi_table_new_column_array (vis_SC, output_name,
"rad", CPL_TYPE_DOUBLE, nwave_sc);
2196 cpl_array ** phaseref = cpl_table_get_data_array (vis_SC, output_name);
2199 gravi_table_new_column_array (vis_SC, coeff_name, NULL, CPL_TYPE_DOUBLE, maxdeg+1);
2200 cpl_array ** phase_coeff = cpl_table_get_data_array (vis_SC, coeff_name);
2202 CPLCHECK_MSG (
"Cannot create column");
2205 for (cpl_size base = 0; base < nbase; base++) {
2206 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2207 int nsc = row_sc * nbase + base;
2211 cpl_array * phaseref_ftdit = cpl_array_cast (visData_ftdit[nsc], CPL_TYPE_DOUBLE_COMPLEX);
2215 double mean_phase = carg (cpl_array_get_mean_complex (phaseref_ftdit));
2216 cpl_array_multiply_scalar_complex (phaseref_ftdit, cexp(- I * mean_phase));
2219 cpl_array_arg (phaseref_ftdit);
2221 cpl_array_add_scalar (phaseref_ftdit, mean_phase);
2227 cpl_vector * input = cpl_vector_wrap (nwave_ft, cpl_array_get_data_double (phaseref_ftdit));
2228 cpl_polynomial_fit (fit, wave_ft, NULL, input, NULL, CPL_FALSE, &mindeg, &maxdeg);
2229 cpl_vector_unwrap (input);
2230 cpl_array_delete (phaseref_ftdit);
2233 phase_coeff[nsc] = cpl_array_new (maxdeg+1, CPL_TYPE_DOUBLE);
2234 for (cpl_size d = 0; d < maxdeg+1; d++)
2235 cpl_array_set (phase_coeff[nsc], d, cpl_polynomial_get_coeff (fit, &d));
2238 phaseref[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2239 for (cpl_size w = 0; w < nwave_sc; w++) {
2240 double delta = (cpl_vector_get(wave_sc, w) - lbd0) / delta0;
2241 cpl_array_set (phaseref[nsc], w, cpl_polynomial_eval_1d (fit, delta, NULL));
2244 gravi_array_phase_wrap (phaseref[nsc]);
2245 cpl_array_multiply_scalar (phaseref[nsc], -1.0);
2247 CPLCHECK_MSG (
"Cannot compute the PHASE_REF for SC");
2251 FREE (cpl_vector_delete, wave_sc);
2252 FREE (cpl_matrix_delete, wave_ft);
2253 FREE (cpl_array_delete, wavenumber_ft);
2254 FREE (cpl_polynomial_delete, fit);
2256 gravi_msg_function_exit(1);
2257 return CPL_ERROR_NONE;
2272 cpl_table * disp_table)
2274 gravi_msg_function_start(1);
2275 cpl_ensure_code (flux_SC, CPL_ERROR_NULL_INPUT);
2278 cpl_size nrow = cpl_table_get_nrow (flux_SC) / ntel;
2281 gravi_table_new_column (flux_SC,
"FDDL",
"m", CPL_TYPE_DOUBLE);
2282 CPLCHECK_MSG (
"Cannot create columns");
2285 if (disp_table == NULL) {
2286 gravi_msg_function_exit(1);
2287 return CPL_ERROR_NONE;
2291 cpl_msg_info (cpl_func,
"Load the linearity model coeficients");
2293 double ** lin_fddl_sc = gravi_table_get_data_array_double (disp_table,
"LIN_FDDL_SC");
2294 double ** lin_fddl_ft = gravi_table_get_data_array_double (disp_table,
"LIN_FDDL_FT");
2295 cpl_size disp_order = cpl_table_get_column_depth (disp_table,
"LIN_FDDL_SC");
2296 CPLCHECK_MSG (
"Cannot get linearity model data");
2299 double * ftpos = cpl_table_get_data_double (flux_SC,
"FT_POS");
2300 double * scpos = cpl_table_get_data_double (flux_SC,
"SC_POS");
2301 double * fddl = cpl_table_get_data_double (flux_SC,
"FDDL");
2302 CPLCHECK_MSG (
"Cannot get POS data");
2305 for (cpl_size tel = 0; tel < ntel; tel++) {
2306 for (cpl_size row = 0; row < nrow; row ++) {
2307 cpl_size nsc = row * ntel + tel;
2312 for (
int o = 0; o < disp_order; o++) {
2313 fddl[nsc] += lin_fddl_sc[tel][o] * pow (scpos[nsc], (
double)o) * 0.5e-6;
2314 fddl[nsc] += lin_fddl_ft[tel][o] * pow (ftpos[nsc], (
double)o) * 0.5e-6;
2320 FREE (cpl_free, lin_fddl_sc);
2321 FREE (cpl_free, lin_fddl_ft);
2323 gravi_msg_function_exit(1);
2324 return CPL_ERROR_NONE;
2344 cpl_table * flux_SC,
2345 cpl_table * wave_table,
2346 cpl_table * disp_table,
2347 cpl_propertylist * header,
2348 const cpl_parameterlist * parlist)
2350 gravi_msg_function_start(1);
2351 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2352 cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
2353 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
2355 cpl_size nwave_sc = cpl_table_get_column_depth (vis_SC,
"VISDATA");
2358 if (disp_table == NULL) {
2359 cpl_msg_info (cpl_func,
"Cannot create OPD_DISP, not DISP_MODEL table");
2360 gravi_msg_function_exit(1);
2361 return CPL_ERROR_NONE;
2367 cpl_msg_info (cpl_func,
"Load the dispersion model coeficients");
2370 double ** n_mean = cpl_malloc (4 *
sizeof(
double*));
2371 double ** n_diff = cpl_malloc (4 *
sizeof(
double*));
2372 for (
int t = 0; t < 4; t++) {
2373 n_mean[t] = cpl_calloc (nwave_sc,
sizeof(
double));
2374 n_diff[t] = cpl_calloc (nwave_sc,
sizeof(
double));
2377 if ( cpl_table_has_column (disp_table,
"N_MEAN") &&
2378 cpl_table_has_column (disp_table,
"WAVE0")) {
2381 cpl_msg_info (cpl_func,
"Dispersion model as N/Nmet [wave0/wave-1]");
2382 cpl_size disp_order = cpl_table_get_column_depth (disp_table,
"N_MEAN");
2383 for (
int t = 0; t < 4; t++) {
2384 double lbd0 = cpl_table_get (disp_table,
"WAVE0", t, NULL);
2385 for (
int w = 0; w < nwave_sc; w++ ) {
2386 double lbd = cpl_table_get (wave_table,
"EFF_WAVE", w, NULL);
2387 double xfit = lbd0/lbd - 1.0;
2388 for (
int o = 0; o < disp_order; o++) {
2389 n_mean[t][w] += gravi_table_get_value (disp_table,
"N_MEAN", t, o) * pow (xfit, o);
2390 n_diff[t][w] += gravi_table_get_value (disp_table,
"N_DIFF", t, o) * pow (xfit, o);
2395 FREELOOP (cpl_free, n_mean, 4);
2396 FREELOOP (cpl_free, n_diff, 4);
2397 cpl_msg_error (cpl_func,
"The DISP_MODEL is not recognized... contact the DRS team");
2398 return cpl_error_set_message (cpl_func,CPL_ERROR_ILLEGAL_INPUT,
2399 "The DISP_MODEL is not recognized... contact the DRS team");
2402 CPLCHECK_MSG (
"Cannot compute N_MEAN or N_DIFF");
2404 cpl_size nbase = 6, ntel = 4;
2405 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
2408 gravi_table_new_column_array (vis_SC,
"OPD_DISP",
"m", CPL_TYPE_DOUBLE, nwave_sc);
2409 cpl_array ** opd_disp = cpl_table_get_data_array (vis_SC,
"OPD_DISP");
2411 gravi_table_new_column (vis_SC,
"GDELAY_DISP",
"m", CPL_TYPE_DOUBLE);
2412 double * gd_disp = cpl_table_get_data_double (vis_SC,
"GDELAY_DISP");
2414 gravi_table_new_column_array (vis_SC,
"PHASE_DISP",
"rad", CPL_TYPE_DOUBLE, nwave_sc);
2415 cpl_array ** phase_disp = cpl_table_get_data_array (vis_SC,
"PHASE_DISP");
2417 CPLCHECK_MSG (
"Cannot create columns");
2420 double * opd_met = cpl_table_get_data_double (flux_SC,
"OPD_MET_FC");
2421 double * fddl = cpl_table_get_data_double (flux_SC,
"FDDL");
2423 CPLCHECK_MSG (
"Cannot get data");
2426 double * wavenumber_sc = cpl_malloc (nwave_sc *
sizeof(
double));
2427 for (cpl_size wave = 0; wave < nwave_sc ; wave ++) {
2428 wavenumber_sc[wave] = 1. / cpl_table_get (wave_table,
"EFF_WAVE", wave, NULL);
2434 double * opl_zero_fc = cpl_malloc (4 *
sizeof(
double));
2438 for (t = 0; t < ntel; t++) {
2439 opl_zero_fc[t] = 0.0;
2441 if ( gravi_param_get_bool (parlist,
"gravity.signal.use-met-zero") ) {
2442 cpl_msg_info (cpl_func,
"Metrology zero calculation is enabled!");
2445 for (t = 0; t < ntel; t++) {
2446 sprintf (name,
"ESO OCS MET OPL_ZERO_FC%i", t+1);
2447 if (cpl_propertylist_has (header, name)) {
2448 opl_zero_fc[t] = cpl_propertylist_get_double (header, name)*1e-3;
2449 sprintf (name,
"ESO FDDL MET OFFSET%i", t+1);
2450 opl_zero_fc[t] += cpl_propertylist_get_double (header, name)*1e-3;
2451 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]);
2456 for (t = 0; t < ntel; t++) {
2457 if ( gravi_pfits_has_gdzero (header, t+1) ) {
2459 gravi_pfits_get_gdzero (header, t+1)*1e-3
2460 * (wavenumber_sc[nwave_sc/2-1] - wavenumber_sc[nwave_sc/2+1])
2461 / (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]);
2462 cpl_msg_info (cpl_func,
"Updating metrology zero with QC/PRO MET GD_ZERO_FC%i: %f [mm]", t+1, opl_zero_fc[t]);
2466 for (t = 0; t < ntel; t++) {
2467 if ( gravi_pfits_get_oplzero (header, t+1) ){
2468 opl_zero_fc[t] = gravi_pfits_get_oplzero (header, t+1)*1e-3;
2469 cpl_msg_info (cpl_func,
"Updating metrology zero with QC/PRO MET OPL_ZERO_FC%i: %f [mm]", t+1, opl_zero_fc[t]);
2473 cpl_msg_info (cpl_func,
"Metrology zero calculation is disabled!");
2477 for (cpl_size base = 0; base < nbase; base++) {
2478 for (cpl_size row_sc = 0; row_sc < nrow_sc; row_sc ++) {
2480 cpl_size nsc = row_sc * nbase + base;
2481 cpl_size t0 = GRAVI_BASE_TEL[base][0], t0f = row_sc * ntel + t0;
2482 cpl_size t1 = GRAVI_BASE_TEL[base][1], t1f = row_sc * ntel + t1;
2486 cpl_size wave = nwave_sc / 2;
2487 double s1 = wavenumber_sc[wave-1];
2488 double s2 = wavenumber_sc[wave+1];
2490 (n_mean[t0][wave-1] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][wave-1] * fddl[t0f]) -
2491 (n_mean[t1][wave-1] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][wave-1] * fddl[t1f]);
2493 (n_mean[t0][wave+1] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][wave+1] * fddl[t0f]) -
2494 (n_mean[t1][wave+1] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][wave+1] * fddl[t1f]);
2496 gd_disp[nsc] = (o1*s1 - o2*s2) / (s1-s2);
2500 phase_disp[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2502 for (
int w = 0; w < nwave_sc ; w++) {
2503 cpl_array_set (phase_disp[nsc], w,
2504 ((n_mean[t0][w] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][w] * fddl[t0f]) -
2505 (n_mean[t1][w] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][w] * fddl[t1f]) -
2506 gd_disp[nsc]) * wavenumber_sc[w] * CPL_MATH_2PI);
2509 gravi_array_phase_wrap (phase_disp[nsc]);
2514 opd_disp[nsc] = cpl_array_new (nwave_sc, CPL_TYPE_DOUBLE);
2516 for (
int w = 0; w < nwave_sc ; w++) {
2517 cpl_array_set (opd_disp[nsc], w,
2518 (n_mean[t0][w] * (opd_met[t0f]-opl_zero_fc[t0]) + n_diff[t0][w] * fddl[t0f]) -
2519 (n_mean[t1][w] * (opd_met[t1f]-opl_zero_fc[t1]) + n_diff[t1][w] * fddl[t1f]) );
2526 FREELOOP (cpl_free, n_mean, 4);
2527 FREELOOP (cpl_free, n_diff, 4);
2528 FREE (cpl_free, wavenumber_sc);
2529 FREE (cpl_free, opl_zero_fc);
2531 gravi_msg_function_exit(1);
2532 return CPL_ERROR_NONE;
2550 cpl_table * wave_table,
2551 cpl_propertylist * header,
2552 const cpl_parameterlist * parlist)
2554 gravi_msg_function_start(1);
2555 cpl_ensure_code (vis_SC, CPL_ERROR_NULL_INPUT);
2556 cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
2557 cpl_ensure_code (header, CPL_ERROR_NULL_INPUT);
2560 cpl_size nwave = cpl_table_get_column_depth (vis_SC,
"VISDATA");
2561 cpl_size nrow = cpl_table_get_nrow (vis_SC);
2564 double * ucoord = cpl_table_get_data_double (vis_SC,
"UCOORD");
2565 double * vcoord = cpl_table_get_data_double (vis_SC,
"VCOORD");
2566 cpl_array ** phase_ref = cpl_table_get_data_array (vis_SC,
"PHASE_REF");
2567 double * opd_met_telfc_mcorr = NULL;
2568 double * opd_met_fc_corr = NULL;
2569 CPLCHECK_MSG (
"Cannot get input data");
2572 if ( !cpl_table_has_column (vis_SC,
"OPD_DISP") ) {
2573 cpl_msg_info (cpl_func,
"Cannot compute IMAGING_REF, not column OPD_DISP");
2574 gravi_msg_function_exit(1);
2575 return CPL_ERROR_NONE;
2578 cpl_array ** opd_disp = cpl_table_get_data_array (vis_SC,
"OPD_DISP");
2579 CPLCHECK_MSG (
"Cannot get OPD_DISP data");
2582 double sep_U = gravi_pfits_get_sobj_x (header)*1e-3/3600.0/CPL_MATH_DEG_RAD;
2583 double sep_V = gravi_pfits_get_sobj_y (header)*1e-3/3600.0/CPL_MATH_DEG_RAD;
2586 gravi_table_new_column_array (vis_SC,
"IMAGING_REF",
"rad", CPL_TYPE_DOUBLE, nwave);
2587 cpl_array ** imaging_ref = cpl_table_get_data_array (vis_SC,
"IMAGING_REF");
2588 CPLCHECK_MSG (
"Cannot create column");
2591 const char * imaging_ref_met = gravi_param_get_string (parlist,
"gravity.signal.imaging-ref-met");
2592 cpl_msg_info (cpl_func,
"imaging-ref-met = %s", imaging_ref_met);
2593 if (!strcmp (imaging_ref_met,
"TEL")) {
2594 cpl_msg_info (cpl_func,
"Use telescope metrology for IMAGING_REF computation");
2595 opd_met_telfc_mcorr = cpl_table_get_data_double (vis_SC,
"OPD_MET_TELFC_MCORR");
2596 opd_met_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
2597 }
else if (!strcmp (imaging_ref_met,
"FC_CORR")) {
2598 cpl_msg_info (cpl_func,
"Use corrected fiber coupler metrology for IMAGING_REF computation");
2599 opd_met_fc_corr = cpl_table_get_data_double (vis_SC,
"OPD_MET_FC_CORR");
2600 }
else if (!strcmp (imaging_ref_met,
"FC")) {
2601 cpl_msg_info (cpl_func,
"Use fiber coupler metrology for IMAGING_REF computation");
2603 cpl_msg_error (cpl_func,
"Unknown metrology source for IMAGING_REF calculation!");
2604 cpl_ensure_code (imaging_ref_met, CPL_ERROR_ILLEGAL_INPUT);
2608 for (cpl_size row = 0; row < nrow; row ++) {
2611 imaging_ref[row] = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
2614 double opd_met_corr = (opd_met_telfc_mcorr ? opd_met_telfc_mcorr[row] : 0.0)
2615 + (opd_met_fc_corr ? opd_met_fc_corr[row] : 0.0);
2618 for (
int w = 0; w < nwave; w++) {
2620 double wavelength = cpl_table_get (wave_table,
"EFF_WAVE", w, NULL);
2624 cpl_array_set (imaging_ref[row], w,
2625 cpl_array_get (phase_ref[row], w, NULL)
2626 - cpl_array_get (opd_disp[row], w, NULL) * CPL_MATH_2PI / wavelength
2627 + (ucoord[row] * sep_U + vcoord[row] * sep_V) * CPL_MATH_2PI / wavelength
2628 - opd_met_corr * CPL_MATH_2PI / wavelength);
2629 CPLCHECK_MSG (
"Cannot compute the imaging phase");
2633 gravi_array_phase_wrap (imaging_ref[row]);
2636 gravi_msg_function_exit(1);
2637 return CPL_ERROR_NONE;
2657 gravi_data * disp_data,
2658 const cpl_parameterlist * parlist)
2660 gravi_msg_function_start(1);
2661 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
2663 int nbase = 6, ntel = 4;
2666 cpl_propertylist * p2vmred_header = gravi_data_get_header (p2vmred_data);
2667 int npol_ft = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_FT);
2668 int npol_sc = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_SC);
2670 CPLCHECK_MSG (
"Cannot get header");
2673 if ( npol_sc != npol_ft ) {
2674 gravi_msg_warning (
"FIXME",
"Not sure this function works with npol_sc != npol_ft");
2679 double dit_sc = gravi_pfits_get_dit_sc (p2vmred_header) * 1e6;
2680 cpl_msg_info (cpl_func,
"dit_sc = %g [us]", dit_sc);
2686 for (
int pol = 0; pol < CPL_MAX(npol_sc,npol_ft); pol++) {
2689 cpl_msg_info (cpl_func,
"Start polarisation %d over %d",pol+1, CPL_MAX(npol_sc,npol_ft));
2690 cpl_msg_info(cpl_func,
"Insname FT : %s, pol %d npol %d",
2691 GRAVI_INSNAME(GRAVI_FT,pol,npol_ft), pol, npol_ft);
2692 cpl_msg_info(cpl_func,
"Insname SC : %s, pol %d npol %d",
2693 GRAVI_INSNAME(GRAVI_SC,pol,npol_sc), pol, npol_sc);
2696 cpl_table * vis_FT = gravi_data_get_oi_vis (p2vmred_data, GRAVI_FT, pol, npol_ft);
2697 cpl_table * flux_FT = gravi_data_get_oi_flux (p2vmred_data, GRAVI_FT, pol, npol_ft);
2698 CPLCHECK_MSG (
"Cannot get the FT tables");
2701 cpl_table * vis_SC = gravi_data_get_oi_vis (p2vmred_data, GRAVI_SC, pol, npol_sc);
2702 cpl_table * flux_SC = gravi_data_get_oi_flux (p2vmred_data, GRAVI_SC, pol, npol_sc);
2703 CPLCHECK_MSG (
"Cannot get the SC tables");
2708 CPLCHECK_MSG (
"Cannot get the VIS_MET and FDDL tables");
2711 cpl_table * oi_wavelengthft = gravi_data_get_oi_wave (p2vmred_data, GRAVI_FT, pol, npol_ft);
2712 cpl_table * oi_wavelengthsc = gravi_data_get_oi_wave (p2vmred_data, GRAVI_SC, pol, npol_sc);
2713 CPLCHECK_MSG (
"Cannot get the OI_WAVELENGTH tables");
2724 CPLCHECK_MSG (
"Cannot sync vis_SC");
2731 CPLCHECK_MSG (
"Cannot sync flux_SC");
2737 cpl_array ** flux_ft = cpl_table_get_data_array (flux_FT,
"FLUX");
2738 cpl_size nwave_ft = cpl_table_get_column_depth (vis_FT,
"VISDATA");
2739 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
2741 gravi_table_new_column (flux_FT,
"TOTALFLUX",
"e", CPL_TYPE_DOUBLE);
2742 double * total_flux_ft = cpl_table_get_data_double (flux_FT,
"TOTALFLUX");
2744 CPLCHECK_MSG (
"Cannot create columns");
2747 for (cpl_size row = 0; row < nrow_ft * ntel; row ++) {
2748 total_flux_ft[row] = cpl_array_get_mean (flux_ft[row]) * nwave_ft;
2753 gravi_table_smooth_column (flux_FT,
"TOTALFLUX",
"TOTALFLUX", 10, ntel);
2755 CPLCHECK_MSG(
"Cannot compute TOTALFLUX for FT");
2765 CPLCHECK_MSG (
"Cannot create signals for VIS_FT");
2774 CPLCHECK_MSG (
"Cannot create the signal for FLUX_SC");
2784 vis_FT, oi_wavelengthft);
2786 CPLCHECK_MSG (
"Cannot create signals for VIS_SC");
2792 for (
int base = 0; base < nbase; base++) {
2795 sprintf (qc_name,
"ESO QC VFACTOR%s_P%d AVG", GRAVI_BASE_NAME[base], pol+1);
2796 double vmean = gravi_table_get_column_mean (vis_SC,
"V_FACTOR_WL", base, nbase);
2797 cpl_propertylist_update_double (p2vmred_header, qc_name, vmean);
2798 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean v-factor");
2800 sprintf (qc_name,
"ESO QC PFACTOR%s_P%d AVG", GRAVI_BASE_NAME[base], pol+1);
2801 double pmean = gravi_table_get_column_mean (vis_SC,
"P_FACTOR", base, nbase);
2802 cpl_propertylist_update_double (p2vmred_header, qc_name, pmean);
2803 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"mean p-factor");
2831 CPLCHECK_MSG (
"Cannot compute the OPD_DISP");
2840 gravi_table_compute_group_delay (vis_FT,
"VISDATA",
"GDELAY", oi_wavelengthft);
2843 gravi_table_compute_group_delay (vis_SC,
"VISDATA",
"GDELAY", oi_wavelengthsc);
2846 gravi_table_compute_group_delay (vis_SC,
"VISDATA_FT",
"GDELAY_FT", oi_wavelengthft);
2848 CPLCHECK_MSG (
"Cannot compute the GDELAYs");
2859 CPLCHECK_MSG (
"Cannot compute the PHASE_REF");
2863 gravi_msg_function_exit(1);
2864 return CPL_ERROR_NONE;
2890 const cpl_parameterlist * parlist)
2892 gravi_msg_function_start(1);
2893 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
2898 cpl_propertylist * p2vmred_header = gravi_data_get_header (p2vmred_data);
2899 CPLCHECK_MSG (
"Cannot get header");
2905 cpl_msg_info (cpl_func,
"Cannot compute rejection flags for FT (no FT data)");
2908 cpl_msg_info (cpl_func,
"Compute rejection flags for FT");
2911 double threshold_SNR_ft = gravi_param_get_double (parlist,
"gravity.signal.snr-min-ft");
2912 double threshold_STATE_ft = gravi_param_get_double (parlist,
"gravity.signal.state-min-ft");
2913 double min_GSTATE_ft = gravi_param_get_double (parlist,
"gravity.signal.global-state-min-ft");
2914 double max_GSTATE_ft = gravi_param_get_double (parlist,
"gravity.signal.global-state-max-ft");
2916 cpl_msg_info (cpl_func,
"SNR threshold to define fringe-detection in FT: %g", threshold_SNR_ft);
2917 cpl_msg_info (cpl_func,
"STATE threshold for FT: %g", threshold_STATE_ft);
2918 cpl_msg_info (cpl_func,
"GLOBAL_STATE threshold for FT: %g - %g", min_GSTATE_ft, max_GSTATE_ft);
2921 int npol_ft = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_FT);
2922 for (
int pol = 0; pol < npol_ft; pol++) {
2925 cpl_table * vis_FT = gravi_data_get_oi_vis (p2vmred_data, GRAVI_FT, pol, npol_ft);
2926 double * snr = cpl_table_get_data_double (vis_FT,
"SNR_BOOT");
2927 int * state = cpl_table_get_data_int (vis_FT,
"STATE");
2928 int * gstate = cpl_table_get_data_int (vis_FT,
"OPDC_STATE");
2930 gravi_table_new_column (vis_FT,
"REJECTION_FLAG", NULL, CPL_TYPE_INT);
2931 int * reject_flag_ft = cpl_table_get_data_int (vis_FT,
"REJECTION_FLAG");
2933 CPLCHECK_MSG (
"Cannot create columns");
2936 cpl_size nrow_ft = cpl_table_get_nrow (vis_FT) / nbase;
2937 for (cpl_size row = 0; row < nrow_ft * nbase; row ++) {
2941 if ( snr[row] < threshold_SNR_ft )
2942 gravi_bit_set (reject_flag_ft[row], snr_bit);
2944 gravi_bit_clear (reject_flag_ft[row], snr_bit);
2948 if ( state[row] < threshold_STATE_ft ||
2949 gstate[row] < min_GSTATE_ft ||
2950 gstate[row] > max_GSTATE_ft )
2951 gravi_bit_set (reject_flag_ft[row], state_bit);
2953 gravi_bit_clear (reject_flag_ft[row], state_bit);
2965 cpl_msg_info (cpl_func,
"Cannot compute tracking ratio for SC (no FT or SC data)");
2968 cpl_msg_info (cpl_func,
"Compute tracking ratio for each SC DIT");
2971 int npol_sc = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_SC);
2972 int npol_ft = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_FT);
2973 for (
int pol = 0; pol < npol_sc; pol++) {
2975 cpl_table * vis_SC = gravi_data_get_oi_vis (p2vmred_data, GRAVI_SC, pol, npol_sc);
2976 cpl_table * vis_FT = gravi_data_get_oi_vis (p2vmred_data, GRAVI_FT, CPL_MIN(pol,npol_ft), npol_ft);
2978 CPLCHECK_MSG (
"Cannot compute lockratio_sc");
2986 cpl_msg_info (cpl_func,
"Cannot compute rejection flags for SC (no SC data)");
2989 cpl_msg_info (cpl_func,
"Compute rejection flags for SC");
2992 double minlockratio_sc = gravi_param_get_double (parlist,
"gravity.signal.tracking-min-sc");
2993 double minvfactor_sc = gravi_param_get_double (parlist,
"gravity.signal.vfactor-min-sc");
2995 cpl_msg_info (cpl_func,
"Fringe-detection ratio to discard frame on SC: %g", minlockratio_sc);
2996 cpl_msg_info (cpl_func,
"vFactor threshold to discard frame on SC: %g", minvfactor_sc);
2999 int npol_sc = gravi_pfits_get_pola_num (p2vmred_header, GRAVI_SC);
3000 for (
int pol = 0; pol < npol_sc; pol++) {
3003 cpl_table * vis_SC = gravi_data_get_oi_vis (p2vmred_data, GRAVI_SC, pol, npol_sc);
3004 double * vFactor_wl = cpl_table_get_data_double (vis_SC,
"V_FACTOR_WL");
3005 double * fringedet_ftdit = cpl_table_get_data_double (vis_SC,
"FRINGEDET_RATIO");
3007 gravi_table_new_column (vis_SC,
"REJECTION_FLAG", NULL, CPL_TYPE_INT);
3008 int * reject_flag_sc = cpl_table_get_data_int (vis_SC,
"REJECTION_FLAG");
3010 CPLCHECK_MSG (
"Cannot create columns");
3013 cpl_size nrow_sc = cpl_table_get_nrow (vis_SC) / nbase;
3014 for (cpl_size row_sc = 0; row_sc < nrow_sc * nbase; row_sc++) {
3018 if ( fringedet_ftdit[row_sc] < minlockratio_sc )
3019 gravi_bit_set (reject_flag_sc[row_sc], lock_bit);
3021 gravi_bit_clear (reject_flag_sc[row_sc], lock_bit);
3024 int vfactor_bit = 1;
3025 if ( vFactor_wl[row_sc] < minvfactor_sc )
3026 gravi_bit_set (reject_flag_sc[row_sc], vfactor_bit);
3028 gravi_bit_clear (reject_flag_sc[row_sc], vfactor_bit);
3034 gravi_msg_function_exit(1);
3035 return CPL_ERROR_NONE;
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_f1f2_ft(cpl_table *vis_FT, cpl_table *flux_FT)
Compute the photometric normalisation for the FT.
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.
int gravi_data_has_extension(gravi_data *raw_calib, const char *ext_name)
Check if data has extension with given EXTNAME.
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_compute_mean_phasor(cpl_table *oi_vis, const char *name_vis, const char *name_err, const char *name_pha, const char *name_var)
Compute real-time mean phasor of a VISDATA by averaging all spectral elements.
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_compute_snr(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Compute real-time SNR and Group-Delay of the observation.
cpl_error_code gravi_vis_create_phaseref_sc(cpl_table *vis_SC, cpl_table *wave_table_sc, cpl_table *wave_table_ft)
Compute the reference phase for each SC DIT.
cpl_error_code gravi_vis_compute_interspectre(cpl_table *oi_vis, const char *name_vis, const char *name_is)
Compute the real-time interspectra.
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_create_phaseref_ft(cpl_table *vis_FT)
Compute the self-reference phase for each FT DIT.
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_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_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_vis_create_pfactor_sc(cpl_table *vis_SC, cpl_table *flux_FT)
Compute the PFACTOR for the SC.
cpl_table * gravi_data_get_table(gravi_data *self, const char *extname)
Return a pointer on a table extension by its EXTNAME.
cpl_error_code gravi_vis_create_met_sc(cpl_table *vis_SC, cpl_table *vis_MET)
Compute the averaged MET signal for each SC DIT per base.
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_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_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_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_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_compute_rejection(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Create rejection flags P2VMREDUCED file.
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_flux_create_fddlpos_sc(cpl_table *flux_SC, cpl_table *fddl_table)
Compute the averaged FDDL signal for each SC DIT per beam.
int gravi_data_has_type(gravi_data *self, const char *type)
Return the number of ext whose EXTNAME and INSNAME match 'type'.
cpl_error_code gravi_array_add_phasors(cpl_array *input, cpl_array *add, cpl_array *sub)
Add a pair of COMPLEX arrays to a COMPLEX array, in-place: input = input + add*conj(sub) ...
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_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_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_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_f1f2_sc(cpl_table *vis_SC, cpl_table *flux_SC)
Compute the photometric normalisation for the SC.