95 cpl_ensure (
header, CPL_ERROR_NULL_INPUT, NULL);
96 cpl_ensure (mode, CPL_ERROR_NULL_INPUT, NULL);
100 if (!(strcmp (mode,
"gravi_dual"))) n_target = 2;
101 else if (!(strcmp (mode,
"gravi_single"))) n_target = 1;
102 else {cpl_error_set_message(cpl_func, CPL_ERROR_NULL_INPUT,
"Invalid mode (bug)");
return NULL;}
105 cpl_table * oi_target = cpl_table_new (n_target);
108 cpl_table_new_column (oi_target,
"TARGET_ID", CPL_TYPE_INT);
109 cpl_table_set_column_savetype (oi_target,
"TARGET_ID", CPL_TYPE_SHORT);
111 cpl_table_new_column (oi_target,
"TARGET", CPL_TYPE_STRING);
113 cpl_table_new_column (oi_target,
"EQUINOX", CPL_TYPE_FLOAT);
114 cpl_table_set_column_unit (oi_target,
"EQUINOX",
"yr");
116 cpl_table_new_column (oi_target,
"RAEP0", CPL_TYPE_DOUBLE);
117 cpl_table_set_column_unit (oi_target,
"RAEP0",
"deg");
118 cpl_table_new_column (oi_target,
"DECEP0", CPL_TYPE_DOUBLE);
119 cpl_table_set_column_unit (oi_target,
"DECEP0",
"deg");
121 cpl_table_new_column (oi_target,
"RA_ERR", CPL_TYPE_DOUBLE);
122 cpl_table_set_column_unit (oi_target,
"RA_ERR",
"deg");
123 cpl_table_new_column (oi_target,
"DEC_ERR", CPL_TYPE_DOUBLE);
124 cpl_table_set_column_unit (oi_target,
"DEC_ERR",
"deg");
126 cpl_table_new_column (oi_target,
"SYSVEL", CPL_TYPE_DOUBLE);
127 cpl_table_set_column_unit (oi_target,
"SYSVEL",
"m/s");
128 cpl_table_new_column (oi_target,
"VELTYP", CPL_TYPE_STRING);
129 cpl_table_new_column (oi_target,
"VELDEF", CPL_TYPE_STRING);
131 cpl_table_new_column (oi_target,
"PMRA", CPL_TYPE_DOUBLE);
132 cpl_table_set_column_unit (oi_target,
"PMRA",
"deg/yr");
133 cpl_table_new_column (oi_target,
"PMDEC", CPL_TYPE_DOUBLE);
134 cpl_table_set_column_unit (oi_target,
"PMDEC",
"deg/yr");
135 cpl_table_new_column (oi_target,
"PMRA_ERR", CPL_TYPE_DOUBLE);
136 cpl_table_set_column_unit (oi_target,
"PMRA_ERR",
"deg/yr");
137 cpl_table_new_column (oi_target,
"PMDEC_ERR", CPL_TYPE_DOUBLE);
138 cpl_table_set_column_unit (oi_target,
"PMDEC_ERR",
"deg/yr");
140 cpl_table_new_column (oi_target,
"PARALLAX", CPL_TYPE_FLOAT);
141 cpl_table_set_column_unit (oi_target,
"PARALLAX",
"deg");
142 cpl_table_new_column (oi_target,
"PARA_ERR", CPL_TYPE_FLOAT);
143 cpl_table_set_column_unit (oi_target,
"PARA_ERR",
"deg");
145 cpl_table_new_column (oi_target,
"SPECTYP", CPL_TYPE_STRING);
149 for (
int ti=0; ti<n_target; ti++) {
152 cpl_table_set_int (oi_target,
"TARGET_ID", ti, ti+1);
171 double dec_err = 0.0;
172 raep *= CPL_MATH_DEG_RAD;
173 decp *= CPL_MATH_DEG_RAD;
177 cpl_msg_info (cpl_func,
"Found target '%s' with ra=%fd and dec=%fd at equinox %.2f", name, raep, decp, equinox);
181 const char * veltyp =
"UNKNOWN ";
182 const char * veldef =
"OPTICAL ";
189 double pmra_err = 0.0;
190 double pmdec_err = 0.0;
195 float para_err = 0.0;
199 const char * spectyp =
"UNKNOWN ";
202 cpl_table_set_float (oi_target,
"EQUINOX", ti, equinox);
203 cpl_table_set_double (oi_target,
"RAEP0", ti, raep);
204 cpl_table_set_double (oi_target,
"DECEP0", ti, decp);
205 cpl_table_set_double (oi_target,
"RA_ERR", ti, ra_err);
206 cpl_table_set_double (oi_target,
"DEC_ERR", ti, dec_err);
207 cpl_table_set_double (oi_target,
"SYSVEL", ti, sysvel);
208 cpl_table_set_string (oi_target,
"VELTYP", ti, veltyp);
209 cpl_table_set_string (oi_target,
"VELDEF", ti, veldef);
210 cpl_table_set_double (oi_target,
"PMRA", ti, pmra);
211 cpl_table_set_double (oi_target,
"PMDEC", ti, pmdec);
212 cpl_table_set_double (oi_target,
"PMRA_ERR", ti, pmra_err);
213 cpl_table_set_double (oi_target,
"PMDEC_ERR",ti, pmdec_err);
214 cpl_table_set_float (oi_target,
"PARALLAX", ti, parallax);
215 cpl_table_set_float (oi_target,
"PARA_ERR", ti, para_err);
216 cpl_table_set_string (oi_target,
"SPECTYP", ti, spectyp);
242 cpl_ensure (array_geometry, CPL_ERROR_NULL_INPUT, NULL);
244 cpl_table * oi_array = cpl_table_duplicate (array_geometry);
247 if (cpl_table_has_column (oi_array,
"STA") ) {
248 cpl_table_name_column (oi_array,
"STA",
"STAXYZ");
253 cpl_size n_row = cpl_table_get_nrow (oi_array);
256 if (is_cal)
cpl_msg_info (cpl_func,
"ARRAY_GEOMETRY has %lld rows instead of 4", n_row);
257 else cpl_msg_warning (cpl_func,
"ARRAY_GEOMETRY has %lld rows instead of 4", n_row);
258 cpl_table_set_size (oi_array, 4);
263 for (
int row=n_row; row<4; row++ ) {
265 cpl_table_set_int (oi_array,
"STA_INDEX", row, 100+row);
266 sprintf (name16,
"S%i", row);
267 cpl_table_set_string (oi_array,
"STA_NAME", row, name16);
268 sprintf (name16,
"T%i", row);
269 cpl_table_set_string (oi_array,
"TEL_NAME", row, name16);
270 cpl_table_set_array (oi_array,
"STAXYZ", row, zero_array);
271 cpl_table_set_float (oi_array,
"DIAMETER", row, 0.0);
272 cpl_table_set_int (oi_array,
"MNTSTA", row, 0);
273 if (is_cal)
cpl_msg_info (cpl_func,
"Add STA_INDEX %i with TEL_NAME=%s in OI_ARRAY",100+row,name16);
274 else cpl_msg_warning (cpl_func,
"Add STA_INDEX %i with TEL_NAME=%s in OI_ARRAY",100+row,name16);
278 FREE (cpl_array_delete, zero_array);
283 for (
int row=0; row<cpl_table_get_nrow (oi_array); row++ ) {
284 const char * tel_name = cpl_table_get_string (oi_array,
"TEL_NAME", row);
286 const char * sta_name = cpl_table_get_string (oi_array,
"STA_NAME", row);
331 const char * mode,
const cpl_parameterlist * parlist)
335 cpl_ensure (preproc_data, CPL_ERROR_NULL_INPUT, NULL);
336 cpl_ensure (p2vm_map, CPL_ERROR_NULL_INPUT, NULL);
337 cpl_ensure (mode, CPL_ERROR_NULL_INPUT, NULL);
338 cpl_ensure (parlist, CPL_ERROR_NULL_INPUT, NULL);
341 int ntel = 4, nclo = 4 ;
345 clock_t timer_ft = 0, timer_sc = 0, timer_other = 0, timer_start = 0;
346 timer_start = clock();
353 if ( (cpl_parameterlist_find_const (parlist,
"gravity.dfs.p2vmred-file") != NULL &&
355 (cpl_parameterlist_find_const (parlist,
"gravity.dfs.astro-file") != NULL &&
357 cpl_msg_info (cpl_func,
"p2vmreduced-file is requested, will add computation time.");
360 cpl_msg_info (cpl_func,
"p2vmreduced-file is not requested.");
373 cpl_propertylist_append (p2vmred_header, preproc_header);
387 CPLCHECK_NUL(
"Cannot convert ARRAY_GEOMETRY into OI_ARRAY");
391 oiarray_plist = cpl_propertylist_duplicate (oiarray_plist);
399 cpl_msg_info (cpl_func,
"Update the ARRAYX/Y/Z keywords as double precision");
410 for (
int lab = 0 ; lab < 4 ; lab++) {
413 cpl_table_unselect_all (optical_train_table);
414 if ( cpl_table_or_selected_int (optical_train_table,
"VALUE2", CPL_EQUAL_TO, labInput[lab]) == 1 ) {
420 int nrow = cpl_table_get_nrow (optical_train_table);
421 cpl_table_set_size (optical_train_table, nrow + 1);
422 cpl_table_set_int (optical_train_table,
"VALUE2", nrow, labInput[lab]);
423 sprintf (name16,
"T%i", nrow);
424 cpl_table_set_string (optical_train_table,
"TEL_NAME", nrow, name16);
425 if (is_cal)
cpl_msg_info (cpl_func,
"Add LABINPUT %i with TEL_NAME=%s in OPTICAL_TRAIN",labInput[lab],name16);
426 else cpl_msg_warning (cpl_func,
"Add LABINPUT %i with TEL_NAME=%s in OPTICAL_TRAIN",labInput[lab],name16);
441 for (
int type_data = 0; type_data < 2; type_data ++){
457 cpl_ensure (cpl_table_get_nrow (detector_table) ==
458 cpl_table_get_nrow (detector_p2vm),
459 CPL_ERROR_ILLEGAL_INPUT, NULL);
462 cpl_size nrow = cpl_table_get_nrow (spectrum_table);
463 int n_detregion = cpl_table_get_nrow (detector_table);
467 int npol = (n_detregion>24 ? 2 : 1);
478 for (
int p = 0; p < npol; p ++) {
480 cpl_msg_info(cpl_func,
"Visibility extraction: polarisation %d over %d for %s",
481 p+1, npol, type_data==
GRAVI_FT?
"FT":
"SC");
483 timer_other += (clock() - timer_start);
484 timer_start = clock();
491 cpl_msg_info(cpl_func,
"Create and set output tables...");
499 cpl_table_erase_column (oi_vis,
"VISAMP");
500 cpl_table_erase_column (oi_vis,
"VISPHI");
501 cpl_table_erase_column (oi_vis,
"VISAMPERR");
502 cpl_table_erase_column (oi_vis,
"VISPHIERR");
503 cpl_table_erase_column (oi_vis,
"RVIS");
504 cpl_table_erase_column (oi_vis,
"RVISERR");
505 cpl_table_erase_column (oi_vis,
"IVIS");
506 cpl_table_erase_column (oi_vis,
"IVISERR");
508 cpl_table_erase_column (oi_t3,
"T3AMP");
509 cpl_table_erase_column (oi_t3,
"T3PHI");
510 cpl_table_erase_column (oi_t3,
"T3AMPERR");
511 cpl_table_erase_column (oi_t3,
"T3PHIERR");
512 cpl_table_erase_column (oi_t3,
"U1COORD");
513 cpl_table_erase_column (oi_t3,
"V1COORD");
514 cpl_table_erase_column (oi_t3,
"U2COORD");
515 cpl_table_erase_column (oi_t3,
"V2COORD");
518 cpl_propertylist * oi_plist = cpl_propertylist_new ();
519 cpl_propertylist_copy_property (oi_plist, oiarray_plist,
"ARRNAME");
520 cpl_propertylist_copy_property (oi_plist, preproc_header,
"DATE-OBS");
521 cpl_propertylist_update_int (oi_plist,
"OI_REVN", 1);
522 cpl_propertylist_update_int (oi_plist,
"NWAVE", nwave);
523 cpl_propertylist_update_string (oi_plist,
"INSNAME",
GRAVI_INSNAME(type_data,p,npol));
524 cpl_propertylist_update_int (oi_plist,
"EXTVER",
GRAVI_EXTVER(type_data,p,npol));
529 oi_plist = cpl_propertylist_duplicate (oi_plist);
534 oi_plist = cpl_propertylist_duplicate (oi_plist);
539 cpl_msg_info(cpl_func,
"Total time to create tables: %.4f s", (
double)(clock()-timer_start)/(
double)CLOCKS_PER_SEC );
546 timer_other += (clock() - timer_start);
547 timer_start = clock();
549 cpl_msg_info (cpl_func,
"Compute the invers V2P2M -> P2VM matrix...");
552 cpl_matrix ** p2vm = cpl_calloc (nwave,
sizeof (cpl_matrix*));
553 cpl_matrix ** v2pm = cpl_calloc (nwave,
sizeof (cpl_matrix*));
558 for (
int detreg = 0; detreg < n_detregion; detreg++) {
560 all_region[n_region] = detreg;
565 const cpl_array ** trans, ** coh, ** phase;
566 trans = cpl_table_get_data_array_const (p2vm_table,
TRANSMISSION);
567 coh = cpl_table_get_data_array_const (p2vm_table,
COHERENCE);
568 phase = cpl_table_get_data_array_const (p2vm_table,
PHASE);
572 double ** pP2VM = cpl_malloc (
sizeof(
double*) * nwave);
573 double ** pV2PM = cpl_malloc (
sizeof(
double*) * nwave);
576 for (cpl_size wave = 0; wave < nwave; wave ++) {
579 v2pm[wave] = cpl_matrix_new (n_region, 16);
582 for (cpl_size region = 0; region < n_region; region++) {
583 int detregion = all_region[region];
587 for (
int i = 0; i < 4; i++){
588 cpl_matrix_set (v2pm[wave], region, i,
589 cpl_array_get (trans[detregion], wave + i * nwave, &nv));
595 for (
int i = 0; i < 6; i++) {
596 cpl_matrix_set (v2pm[wave], region, i + 4,
597 cpl_array_get (coh[detregion], wave + i * nwave, &nv) *
598 cos(cpl_array_get (phase[detregion], wave + i * nwave, &nv)));
599 cpl_matrix_set (v2pm[wave], region, i + 10,
600 cpl_array_get (coh[detregion], wave + i * nwave, &nv) *
601 sin(cpl_array_get (phase[detregion], wave + i * nwave, &nv)));
607 cpl_matrix_multiply_scalar (v2pm[wave], 1./n_region);
614 pP2VM[wave] = cpl_matrix_get_data (p2vm[wave]);
615 pV2PM[wave] = cpl_matrix_get_data (v2pm[wave]);
621 cpl_msg_info (cpl_func,
"Total time to invers matrix: %.4f s", (
double)(clock()-timer_start)/(
double)CLOCKS_PER_SEC );
626 cpl_msg_info (cpl_func,
"Fill the STA_INDEX and TIME.");
627 timer_other += (clock() - timer_start);
628 timer_start = clock();
632 cpl_array * times = cpl_array_wrap_int (cpl_table_get_data_int (spectrum_table,
"TIME"),
633 cpl_table_get_nrow (spectrum_table));
634 cpl_ensure (times, CPL_ERROR_ILLEGAL_INPUT, NULL);
638 cpl_array * mjds = cpl_array_cast (times, CPL_TYPE_DOUBLE);
639 cpl_array_divide_scalar (mjds, 86400.E6);
641 cpl_array_add_scalar (mjds, mjd0);
646 cpl_table_set_column_unit (oi_vis,
"TIME",
"us");
647 cpl_table_set_column_unit (oi_t3,
"TIME",
"us");
648 cpl_table_set_column_unit (oi_flux,
"TIME",
"us");
651 cpl_array * sta_index = cpl_array_new (2, CPL_TYPE_INT);
657 cpl_array_set_int (sta_index, 0, sta0);
658 cpl_array_set_int (sta_index, 1, sta1);
662 for (
int row=0; row < nrow; ++row) {
664 cpl_table_set_array (oi_vis,
"STA_INDEX", idx, sta_index);
665 cpl_table_set (oi_vis,
"TIME", idx, cpl_array_get (times, row, &nv));
666 cpl_table_set (oi_vis,
"MJD", idx, cpl_array_get (mjds, row, &nv));
670 cpl_array_delete (sta_index);
671 CPLCHECK_NUL (
"Cannot fill sta_index or time in OI_VIS");
674 sta_index = cpl_array_new (3, CPL_TYPE_INT);
676 for (
int closure=0; closure < nclo; ++closure) {
681 cpl_array_set_int (sta_index, 0, sta0);
682 cpl_array_set_int (sta_index, 1, sta1);
683 cpl_array_set_int (sta_index, 2, sta2);
687 for (
int row=0; row < nrow; ++row) {
688 int idx = row * nclo + closure;
689 cpl_table_set_array (oi_t3,
"STA_INDEX", idx, sta_index);
690 cpl_table_set (oi_t3,
"TIME", idx, cpl_array_get (times, row, &nv));
691 cpl_table_set (oi_t3,
"MJD", idx, cpl_array_get (mjds, row, &nv));
695 cpl_array_delete (sta_index);
696 CPLCHECK_NUL (
"Cannot fill sta_index or time in OI_T3");
699 for (
int tel = 0; tel <
ntel; tel++){
701 for (
int row=0; row < nrow; ++row) {
702 int idx = row *
ntel + tel;
703 cpl_table_set_int (oi_flux,
"STA_INDEX", idx, sta0);
704 cpl_table_set (oi_flux,
"TIME", idx, cpl_array_get (times, row, &nv));
705 cpl_table_set (oi_flux,
"MJD", idx, cpl_array_get (mjds, row, &nv));
709 cpl_array_delete (mjds);
710 cpl_array_unwrap (times);
711 CPLCHECK_NUL (
"Cannot fill sta_index or time in OI_FLUX");
715 cpl_table_fill_column_window (oi_vis,
"INT_TIME", 0, nrow *
GRAVI_NBASE, exptime);
716 cpl_table_fill_column_window (oi_t3,
"INT_TIME", 0, nrow * nclo, exptime);
717 cpl_table_fill_column_window (oi_flux,
"INT_TIME", 0, nrow *
ntel, exptime);
722 int target_id = (!strcmp(mode,
"gravi_dual") && (type_data==
GRAVI_SC) )?2:1;
723 cpl_table_fill_column_window_int (oi_vis,
"TARGET_ID", 0, nrow *
GRAVI_NBASE, target_id);
724 cpl_table_fill_column_window_int (oi_t3,
"TARGET_ID", 0, nrow * nclo, target_id);
725 cpl_table_fill_column_window_int (oi_flux,
"TARGET_ID", 0, nrow *
ntel, target_id);
728 cpl_msg_info (cpl_func,
"Total time to fill STA_INDEX and TIME: %.4f s",
729 (
double)(clock()-timer_start)/(
double)CLOCKS_PER_SEC );
738 timer_other += (clock() - timer_start);
739 timer_start = clock();
742 double ** pReg = cpl_malloc (n_region *
sizeof(
double*));
743 double ** pErr = cpl_malloc (n_region *
sizeof(
double*));
744 cpl_array *** pRegArr = cpl_malloc (n_region *
sizeof(cpl_array**));
745 cpl_array *** pErrArr = cpl_malloc (n_region *
sizeof(cpl_array**));
746 for (
int reg = 0; reg < n_region; reg++) {
747 pRegArr[reg] = cpl_table_get_data_array (spectrum_table,
GRAVI_DATA[all_region[reg]]);
748 pErrArr[reg] = cpl_table_get_data_array (spectrum_table,
GRAVI_DATAERR[all_region[reg]]);
753 cpl_array** tFlux = cpl_table_get_data_array (oi_flux,
"FLUX");
754 cpl_array** tFluxErr = cpl_table_get_data_array (oi_flux,
"FLUXERR");
755 cpl_array** tVis = cpl_table_get_data_array (oi_vis,
"VISDATA");
756 cpl_array** tVisErr = cpl_table_get_data_array (oi_vis,
"VISERR");
761 cpl_array** tChi2 = cpl_table_get_data_array (oi_flux,
"CHI2");
762 int ndof = n_region - 16;
765 double* pOut = cpl_malloc (16 * nwave *
sizeof(
double));
766 double* pOutVar = cpl_malloc (16 * nwave *
sizeof(
double));
767 double* pChi2 = cpl_malloc (nwave *
sizeof(
double));
770 double full_flux_reg = 0.0, full_flux_tel = 0.0;
775 for (cpl_size row = 0; row < nrow; row ++) {
778 for (
int reg = 0; reg < n_region; reg++) {
779 pReg[reg] = cpl_array_get_data_double (pRegArr[reg][row]);
780 pErr[reg] = cpl_array_get_data_double (pErrArr[reg][row]);
784 for (cpl_size wave = 0 ; wave < nwave ; wave++ ) {
787 for (
int reg = 0; reg < n_region; reg++) full_flux_reg += pReg[reg][wave];
793 for (
int out=0; out<16; out++) {
794 pOut[out*nwave+wave] = 0.0;
795 pOutVar[out*nwave+wave] = 0.0;
796 for (cpl_size reg = 0; reg < n_region; reg++) {
797 pOut[out*nwave+wave] += pReg[reg][wave] * pP2VM[wave][out*n_region+reg];
798 pOutVar[out*nwave+wave] +=
gravi_pow2 (pErr[reg][wave] * pP2VM[wave][out*n_region+reg]);
805 for (cpl_size reg = 0; reg < n_region; reg++) {
807 for (
int out=0; out<16; out++)
808 value += pV2PM[wave][reg*16+out] * pOut[out*nwave+wave];
809 pChi2[wave] +=
gravi_pow2 ((value-pReg[reg][wave]) / pErr[reg][wave]) / ndof;
813 for (
int tel = 0; tel < 4; tel++) full_flux_tel += pOut[tel*nwave+wave];
819 for (
int tel = 0; tel <
ntel; tel++){
820 double * data = cpl_malloc (nwave *
sizeof(
double));
821 for (cpl_size wave = 0 ; wave < nwave ; wave++ )
822 data[wave] = pChi2[wave];
823 tChi2[row*
ntel+tel] = cpl_array_wrap_double (data, nwave);
827 for (
int tel = 0; tel <
ntel; tel++){
828 double * data = cpl_malloc (nwave *
sizeof(
double));
829 for (cpl_size wave = 0 ; wave < nwave ; wave++ )
830 data[wave] = pOut[tel*nwave+wave];
831 tFlux[row*
ntel+tel] = cpl_array_wrap_double (data, nwave);
835 for (
int tel = 0; tel <
ntel; tel++){
836 double * data = cpl_malloc (nwave *
sizeof(
double));
837 for (cpl_size wave = 0 ; wave < nwave ; wave++ )
838 data[wave] = sqrt (pOutVar[tel*nwave+wave]);
839 tFluxErr[row*
ntel+tel] = cpl_array_wrap_double (data, nwave);
844 double complex * data = cpl_malloc (nwave *
sizeof(
double complex));
845 for (cpl_size wave = 0 ; wave < nwave ; wave++ )
846 data[wave] = (
double complex)( pOut[(base+4)*nwave+wave] + 1.*I * pOut[(base+10)*nwave+wave] );
847 tVis[row*
GRAVI_NBASE+base] = cpl_array_wrap_double_complex (data, nwave);
852 double complex * data = cpl_malloc (nwave *
sizeof(
double complex));
853 for (cpl_size wave = 0 ; wave < nwave ; wave++ )
854 data[wave] = (
double complex)( sqrt(pOutVar[(base+4)*nwave+wave]) + 1.*I * sqrt(pOutVar[(base+10)*nwave+wave]) );
855 tVisErr[row*
GRAVI_NBASE+base] = cpl_array_wrap_double_complex (data, nwave);
860 for (cpl_size reg = 0; reg < n_region; reg++) {
861 FREE (cpl_array_delete, pRegArr[reg][row]);
862 FREE (cpl_array_delete, pErrArr[reg][row]);
870 FREE (cpl_free, pReg);
871 FREE (cpl_free, pErr);
872 FREE (cpl_free, pRegArr);
873 FREE (cpl_free, pErrArr);
874 FREE (cpl_free, pOut);
875 FREE (cpl_free, pOutVar);
876 FREE (cpl_free, pChi2);
877 FREE (cpl_free, pP2VM);
878 FREE (cpl_free, pV2PM);
879 FREELOOP (cpl_matrix_delete, p2vm, nwave);
880 FREELOOP (cpl_matrix_delete, v2pm, nwave);
883 cpl_msg_info (cpl_func,
"Total flux in TELs: %.2f [e], in REGIONs:%.2f [e] (ratio=%.5f)",
884 full_flux_tel,full_flux_reg,full_flux_tel/full_flux_reg);
886 sprintf (qc_name,
"ESO QC TRANS P2VM %s",
GRAVI_TYPE(type_data));
888 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[e/e] at P2VM extraction");
891 cpl_msg_info (cpl_func,
"Total time to apply matrix: %.4f s", (
double)(clock()-timer_start)/(
double)CLOCKS_PER_SEC );
894 timer_ft += (clock() - timer_start);
896 timer_sc += (clock() - timer_start);
897 timer_start = clock();
908 sprintf (qc_name,
"ESO QC MET LAMBDA MEAN");
910 cpl_propertylist_set_comment (p2vmred_header, qc_name,
"[m] Effective mean metrology wavelength");
917 timer_other += (clock() - timer_start);
918 cpl_msg_info (cpl_func,
"Total time for FT: %10.4f s", (
double)(timer_ft)/(
double)CLOCKS_PER_SEC);
919 cpl_msg_info (cpl_func,
"Total time for SC: %10.4f s", (
double)(timer_sc)/(
double)CLOCKS_PER_SEC);
920 cpl_msg_info (cpl_func,
"Total time for OTHER: %7.4f s", (
double)(timer_other)/(
double)CLOCKS_PER_SEC);
945 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
956 cpl_size nrow_ft = cpl_table_get_nrow (oi_vis) /
GRAVI_NBASE;
957 cpl_size nrow_opdc = cpl_table_get_nrow (opdc);
964 cpl_table_new_column (oi_flux,
"STATE", CPL_TYPE_INT);
965 cpl_table_fill_column_window (oi_flux,
"STATE", 0, nrow_ft *
ntel, -1);
968 cpl_table_new_column (oi_vis,
"STATE", CPL_TYPE_INT);
969 cpl_table_fill_column_window (oi_vis,
"STATE", 0, nrow_ft *
GRAVI_NBASE, -1);
972 cpl_table_new_column (oi_vis,
"OPDC_STATE", CPL_TYPE_INT);
973 cpl_table_fill_column_window (oi_vis,
"OPDC_STATE", 0, nrow_ft *
GRAVI_NBASE, -1);
975 if (nrow_opdc < nrow_ft)
976 cpl_msg_warning (cpl_func,
"Missing FT or OPDC data: nrow_ft - nrow_opdc = %lli", nrow_ft-nrow_opdc);
980 int * time_opdc = cpl_table_get_data_int (opdc,
"TIME");
981 double * time_ft = cpl_table_get_data_double (oi_vis,
"TIME");
982 int * global_state_opdc = cpl_table_get_data_int (opdc,
"STATE");
983 int * global_state = cpl_table_get_data_int (oi_vis,
"OPDC_STATE");
987 for (cpl_size row_opdc=0, row_ft=0 ; row_ft<nrow_ft ; row_ft++) {
991 if ( (time_ft[row_ft*
GRAVI_NBASE] < time_opdc[0]) || (time_ft[row_ft*
GRAVI_NBASE] > time_opdc[nrow_opdc-1]) )
continue;
992 while ( time_ft[row_ft*
GRAVI_NBASE] > time_opdc[row_opdc] ) row_opdc ++;
996 global_state[row_ft*
GRAVI_NBASE+base] = global_state_opdc[row_opdc];
1000 if ( cpl_table_has_column (opdc,
"BASELINE_STATE") ) {
1002 int * steps_opdc = cpl_table_get_data_int (opdc,
"STEPS");
1003 int * state_opdc = cpl_table_get_data_int (opdc,
"BASELINE_STATE");
1004 int * base_state = cpl_table_get_data_int (oi_vis,
"STATE");
1005 int * tel_state = cpl_table_get_data_int (oi_flux,
"STATE");
1006 double * base_steps = cpl_table_get_data_double (oi_vis,
"TARGET_PHASE");
1010 for (cpl_size row_opdc=0, row_ft=0 ; row_ft<nrow_ft ; row_ft++) {
1014 if ( (time_ft[row_ft*
GRAVI_NBASE] < time_opdc[0]) || (time_ft[row_ft*
GRAVI_NBASE] > time_opdc[nrow_opdc-1]) )
continue;
1015 while ( time_ft[row_ft*
GRAVI_NBASE] > time_opdc[row_opdc] ) row_opdc ++;
1018 int state_flag = state_opdc[row_opdc];
1022 for (
int tel = 0; tel <
ntel; tel++)
1037 for (cpl_size base = 0; base <
GRAVI_NBASE; base++) {
1051 for (
int tel = 0; tel <
ntel; tel++) {
1052 int pos = (tel<3) ? 4*tel : 4*(tel+1);
1053 tmp_arr[tel] = 15 & ((steps_opdc[row_opdc]) >> pos);
1064 cpl_msg_warning (cpl_func,
"No column BASELINE_STATE in OPDC... old data ?");
1065 cpl_msg_warning (cpl_func,
"The STATE flags are set to 1 (valid) although the information is unavailable");
1066 cpl_table_fill_column_window (oi_flux,
"STATE", 0, nrow_ft *
ntel, 1);
1067 cpl_table_fill_column_window (oi_vis,
"STATE", 0, nrow_ft *
GRAVI_NBASE, 1);
1073 cpl_msg_debug (cpl_func,
"Duplicate the FT tracking state for 2sd polarisation");
1076 cpl_table_duplicate_column (tmp,
"TARGET_PHASE", oi_vis,
"TARGET_PHASE");
1077 cpl_table_duplicate_column (tmp,
"STATE", oi_vis,
"STATE");
1078 cpl_table_duplicate_column (tmp,
"OPDC_STATE", oi_vis,
"OPDC_STATE");
1081 cpl_table_duplicate_column (tmp,
"STATE", oi_flux,
"STATE");
1089 cpl_array **visdata = cpl_table_get_data_array (oi_vis,
"VISDATA");
1090 int * state = cpl_table_get_data_int (oi_vis,
"STATE");
1091 double * target_phase = cpl_table_get_data_double (oi_vis,
"TARGET_PHASE");
1094 double complex * tmp_cpx = cpl_malloc (
sizeof(
double complex)*nrow_ft);
1100 double complex mean_cpx = 0.0 + I * 0.0;
1101 for (cpl_size row_ft=0 ; row_ft<nrow_ft ; row_ft++) {
1103 tmp_cpx[row_ft] = cpl_array_get_mean_complex (visdata[row_ft*
GRAVI_NBASE+base]) * cexp (-I*target_phase[row_ft*
GRAVI_NBASE+base]);
1104 mean_cpx += tmp_cpx[row_ft];
1108 double sum = 0.0000001, sum2 = 0.0;
1109 for (cpl_size row_ft=0 ; row_ft<nrow_ft ; row_ft++) {
1111 sum2 += pow (carg (tmp_cpx[row_ft] * conj(mean_cpx)), 2);
1118 cpl_propertylist_set_comment (
header, qc_name,
"[rad] residuals when tracking");
1122 cpl_propertylist_update_int (
header, qc_name, sum*100.0/nrow_ft);
1123 cpl_propertylist_set_comment (
header, qc_name,
"[%] ratio of time with tracking");
1127 FREE (cpl_free,tmp_cpx);
1137 for (cpl_size row_opdc=0; row_opdc < nrow_opdc; row_opdc++) {
1138 if (tab_state[row_opdc] == 3 || tab_state[row_opdc] == 2) tracking ++;
1141 sprintf (qc_name,
"ESO QC TRACKING_RATIO");
1142 cpl_propertylist_update_int (
header, qc_name, tracking*100/nrow_ft);
1143 cpl_propertylist_set_comment (
header, qc_name,
"[%] ratio of time with full FT tracking");
1149 return CPL_ERROR_NONE;
1169 cpl_ensure_code (data, CPL_ERROR_NULL_INPUT);
1177 cpl_size nrow = cpl_table_get_nrow (opdc);
1182 cpl_size max_delay = 200;
1185 int *time = cpl_table_get_data_int (opdc,
"TIME");
1197 sprintf (qc_name,
"ESO QC TAU0 OPDC%d%d", t1+1,t2+1);
1199 cpl_propertylist_set_comment (
header, qc_name,
"[s] tau0 for variance of 1 rad2");
1203 for (cpl_size delay=0; delay < max_delay; delay++) {
1207 for (cpl_size row=0; row < nrow - max_delay - 1 ; row++) {
1208 float diff = (piezo[row][t1]-piezo[row][t2]) - (piezo[row+delay][t1]-piezo[row+delay][t2]);
1209 strf += (double)(diff * diff);
1211 strf = strf * gain * gain / (nrow - max_delay - 1);
1216 double tau0 = (time[delay+1] - time[1]) * 1e-6;
1217 cpl_msg_info (cpl_func,
"Compute %s = %.2f [ms]", qc_name, 1e3 * tau0);
1227 FREE (cpl_free, piezo);
1230 return CPL_ERROR_NONE;
1248 cpl_ensure_code (data, CPL_ERROR_NULL_INPUT);
1251 double p05 = 0.0, p95 = 0.0;
1255 int ntel = 4, hw = 2;
1258 cpl_msg_info (cpl_func,
"Smooth flux sery over %i samples", hw);
1262 cpl_size nrow = cpl_table_get_nrow (table) /
ntel;
1266 for (
int tel = 0; tel<
ntel; tel++) {
1269 cpl_vector * raw_flux = cpl_vector_new (nrow);
1270 cpl_vector_fill (raw_flux, 0.0);
1272 for (
int pol = 0; pol<npol; pol++) {
1276 const cpl_array ** flux_array = cpl_table_get_data_array_const (table,
"FLUX");
1278 for (
int row = 0; row<nrow; row++)
1279 cpl_vector_set (raw_flux, row, cpl_vector_get (raw_flux, row) +
1280 cpl_array_get_mean (flux_array[tel+row*
ntel]));
1285 flux = cpl_vector_filter_lowpass_create (raw_flux, CPL_LOWPASS_GAUSSIAN, hw);
1286 cpl_vector_delete (raw_flux);
1289 for (
int row = 0; row<nrow; row++)
1290 cpl_vector_set (flux, row, CPL_MAX (0.0, cpl_vector_get (flux, row)));
1293 cpl_vector_sort (flux, CPL_SORT_ASCENDING);
1295 if (cpl_vector_get_max (flux) > 0)
1296 cpl_vector_divide_scalar (flux, cpl_vector_get_max (flux) / 100);
1300 sprintf (qc_value,
"%.0f,%.0f,%.0f,%.0f,%.0f,%.0f,%.0f,%.0f,%.0f",
1301 cpl_vector_get (flux, (cpl_size)(0.1*(nrow-1))),
1302 cpl_vector_get (flux, (cpl_size)(0.2*(nrow-1))),
1303 cpl_vector_get (flux, (cpl_size)(0.3*(nrow-1))),
1304 cpl_vector_get (flux, (cpl_size)(0.4*(nrow-1))),
1305 cpl_vector_get (flux, (cpl_size)(0.5*(nrow-1))),
1306 cpl_vector_get (flux, (cpl_size)(0.6*(nrow-1))),
1307 cpl_vector_get (flux, (cpl_size)(0.7*(nrow-1))),
1308 cpl_vector_get (flux, (cpl_size)(0.8*(nrow-1))),
1309 cpl_vector_get (flux, (cpl_size)(0.9*(nrow-1))));
1313 sprintf (qc_name,
"ESO QC FLUX_FT%d HISTO", tel+1);
1314 cpl_propertylist_update_string (
header, qc_name, qc_value);
1315 cpl_propertylist_set_comment (
header, qc_name,
"decile in percent");
1318 p05 = cpl_vector_get (flux, (cpl_size)(0.05*(nrow-1)));
1319 p95 = cpl_vector_get (flux, (cpl_size)(0.95*(nrow-1)));
1324 sprintf (qc_name,
"ESO QC FLUX_FT%d P05P95", tel+1);
1326 cpl_propertylist_set_comment (
header, qc_name,
"injected flux 5 percentile over 95 percentile");
1329 cpl_vector_delete (flux);
1333 return CPL_ERROR_NONE;
1352 cpl_ensure_code (p2vmred_data, CPL_ERROR_NULL_INPUT);
1361 int * time_opdc = cpl_table_get_data_int (opdc,
"TIME");
1362 if (cpl_table_has_column(opdc,
"OPD") == 1){
1363 cpl_array **opd_data_opdc = cpl_table_get_data_array (opdc,
"OPD");
1364 cpl_size nrow = cpl_table_get_nrow (opdc);
1368 cpl_array * opd_res_array = cpl_array_new(nrow, CPL_TYPE_DOUBLE_COMPLEX);
1369 cpl_array * opd_array = cpl_array_new(nrow, CPL_TYPE_DOUBLE);
1370 cpl_array * snr_array = cpl_array_new(nrow, CPL_TYPE_DOUBLE);
1373 for (
int pol = 0; pol<npol; pol++) {
1375 double * time_oivis = cpl_table_get_data_double (oi_vis,
"TIME");
1376 cpl_array **visdata = cpl_table_get_data_array (oi_vis,
"VISDATA");
1377 cpl_array **visdata_err = cpl_table_get_data_array (oi_vis,
"VISERR");
1379 cpl_size nrow_oivis = cpl_table_get_nrow (oi_vis);
1380 cpl_size nwave= cpl_array_get_size (visdata[0]);
1385 cpl_size row_oivis = base;
1387 for (cpl_size row = 0; row<nrow; row++)
1391 while (( fabs (time_oivis[row_oivis+
GRAVI_NBASE]-time_opdc[row]+0.0011/2*1e6 )< fabs (time_oivis[row_oivis]-time_opdc[row]+0.0011/2*1e6 ) ) && (row_oivis<nrow_oivis-2*
GRAVI_NBASE))
1396 double complex visdata_mean=cpl_array_get_mean_complex (visdata[row_oivis]);
1397 cpl_array_set_double_complex(opd_res_array,row, visdata_mean*cexp ( I * cpl_array_get_float (opd_data_opdc[row], base, &nv) ));
1398 cpl_array_set_double(opd_array, row, carg(visdata_mean));
1404 for (cpl_size wave = 0 ; wave < nwave ; wave ++)
1406 signal+=cabs(cpl_array_get_double_complex (visdata[row_oivis],wave, &nv));
1407 double err=cabs(cpl_array_get_double_complex (visdata_err[row_oivis],wave, &nv));
1410 noise=sqrt(noise+1);
1412 cpl_array_set_double(snr_array,row, signal /noise -1 );
1414 cpl_array_set_double(snr_array,row, 0 );
1417 CPLCHECK_MSG (
"Cannot compute new array for cross-referencing OPD and FT OPD");
1419 double complex opd_mean = cpl_array_get_mean_complex (opd_res_array);
1420 cpl_array_multiply_scalar_complex(opd_res_array,conj(opd_mean));
1423 double a_num=0,b_num=0,c_num=0,d_num=0,a_den=1,b_den=1,c_den=1,d_den=1;
1424 double opd_std_num=0,opd_std_den=1;
1425 double opd, opd_ref, snr_ref;
1427 for (cpl_size row = 0; row<nrow; row++)
1429 opd=cpl_array_get_double(opd_array,row,&nv);
1430 opd_ref=carg(cpl_array_get_double_complex(opd_res_array,row,&nv));
1431 snr_ref=cpl_array_get_double(snr_array,row,&nv);
1433 a_num+=snr_ref*opd_ref*c_o;
1434 a_den+=snr_ref*c_o*c_o;
1436 b_num+=snr_ref*opd_ref*s_o;
1437 b_den+=snr_ref*s_o*s_o;
1439 c_num+=snr_ref*opd_ref*c_o;
1440 c_den+=snr_ref*c_o*c_o;
1442 d_num+=snr_ref*opd_ref*s_o;
1443 d_den+=snr_ref*s_o*s_o;
1446 opd_std_num+=snr_ref*opd_ref*snr_ref*opd_ref;
1447 opd_std_den+=snr_ref;
1451 CPLCHECK_MSG (
"Cannot calculate variance of linearity");
1453 double a=a_num/a_den;
1454 double b=b_num/b_den;
1455 double c=c_num/c_den;
1456 double d=d_num/d_den;
1458 double opd_std=sqrt(0.1+opd_std_num)/opd_std_den;
1461 double res_car=a*a+b*b+c*c+d*d-16*opd_std*opd_std;
1462 if (res_car > 0) res_car=sqrt(res_car);
1466 if (npol == 2) sprintf (qc_name,
"ESO QC LIN_FT P%d_B%d", pol,base+1);
1467 if (npol == 1) sprintf (qc_name,
"ESO QC LIN_FT P%d_B%d", 3,base+1);
1469 cpl_propertylist_set_comment (
header, qc_name,
"FT nonlinearity biases [rad]");
1476 cpl_array_delete (opd_array);
1477 cpl_array_delete (opd_res_array);
1478 cpl_array_delete (snr_array);
1482 cpl_msg_info(cpl_func,
"Cannot found the OPD column (probably old data). Skip computation of QC LIN_FT");
1486 return CPL_ERROR_NONE;
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
#define gravi_data_get_spectrum_data(data, type)
#define gravi_data_get_oi_flux(data, type, pol, npol)
#define gravi_data_get_header(data)
#define gravi_data_get_p2vm_data(data, type)
#define gravi_data_get_imaging_detector(data, type)
#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_INSNAME(type, pol, npol)
#define GRAVI_OI_TARGET_EXT
#define GRAVI_OI_ARRAY_EXT
#define GRAVI_OI_FLUX_EXT
#define gravi_pfits_get_dit(plist, type)
#define GRAVI_EXTVER(type, pol, npol)
#define gravi_pfits_get_sobj_decep(plist)
#define gravi_pfits_get_sobj_raep(plist)
#define GRAVI_SPECTRUM_DATA_EXT(type)
#define GRAVI_ARRAY_GEOMETRY_EXT
#define GRAVI_METROLOGY_EXT
#define GRAVI_OI_WAVELENGTH_EXT
#define GRAVI_OPTICAL_TRAIN_EXT
#define gravi_bit_get(number, pos)
#define gravi_msg_function_exit(flag)
#define FREE(function, variable)
#define CPLCHECK_NUL(msg)
#define gravi_msg_function_start(flag)
#define CPLCHECK_MSG(msg)
#define FREELOOP(function, variable, n)
cpl_error_code gravi_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
cpl_matrix * gravi_matrix_invertSV_create(cpl_matrix *a_in)
Invers a matrix with singular value decomposition.
cpl_array * gravi_array_init_double(long n, double value)
cpl_error_code gravi_table_set_string_fixlen(cpl_table *table, const char *name, int row, const char *value, int len)
Set string in table, ensuring fixed length (right space padding). see cpl_table_set_string.
float ** gravi_table_get_data_array_float(cpl_table *table, const char *name)
cpl_error_code gravi_table_new_column_array(cpl_table *table, const char *name, const char *unit, cpl_type type, cpl_size size)
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
gravi_data * gravi_data_new(int nb_ext)
Create an empty gravi_data.
cpl_error_code gravi_data_add_table(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_table *table)
Add a BINTABLE extension in gravi_data.
int gravi_data_has_extension(gravi_data *raw_calib, const char *ext_name)
Check if data has extension with given EXTNAME.
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_data_copy_ext(gravi_data *output, gravi_data *input, const char *name)
Copy extensions from one data to another.
int gravi_param_get_bool(const cpl_parameterlist *parlist, const char *name)
cpl_error_code gravi_compute_qc_injection(gravi_data *data)
Compute the QC for the injection stability.
cpl_error_code gravi_compute_tau0(gravi_data *data)
Compute the QC TAU0 parameter.
cpl_error_code gravi_compute_opdc_state(gravi_data *p2vmred_data)
Compute the real-time tracking state from OPDC.
gravi_data * gravi_compute_p2vmred(gravi_data *preproc_data, gravi_data *p2vm_map, const char *mode, const cpl_parameterlist *parlist)
Converts preprocessed data into coherent fluxes using the P2VM.
cpl_table * gravi_create_oitarget_table(const cpl_propertylist *header, const char *mode)
Create the OI_TARGET table from the main header.
cpl_table * gravi_create_oiarray_table(const cpl_table *array_geometry, int is_cal)
Create the OI_ARRAY table from the ARRAY_GEOMETRY.
cpl_error_code gravi_compute_qc_ft_opd_estimator(gravi_data *p2vmred_data)
Compute the QC for the FT linearity.
int gravi_pfits_get_pola_num(const cpl_propertylist *plist, int type_data)
const char * gravi_pfits_get_sobj(const cpl_propertylist *plist)
double gravi_convert_to_mjd(const char *start)
double gravi_pfits_get_pmra(const cpl_propertylist *plist)
cpl_error_code gravi_pfits_ensure_double(cpl_propertylist *self, const char *name)
double gravi_pfits_get_plx(const cpl_propertylist *plist)
int gravi_pfits_is_calib(const cpl_propertylist *plist)
double gravi_pfits_get_robj_raep(const cpl_propertylist *plist)
const char * gravi_pfits_get_robj(const cpl_propertylist *plist)
double gravi_pfits_get_met_wavelength_mean(const cpl_propertylist *plist, cpl_table *met_table)
const char * gravi_pfits_get_start_prcacq(const cpl_propertylist *plist)
double gravi_pfits_get_pmdec(const cpl_propertylist *plist)
double gravi_pfits_get_double_silentdefault(const cpl_propertylist *plist, const char *name, double def)
double gravi_pfits_get_robj_decep(const cpl_propertylist *plist)
cpl_error_code gravi_dump_the_boss(double ra, double dec)
int gravi_region_get_pol(cpl_table *imaging_detector, int region)
Return the polarisation id of a region.
cpl_size gravi_spectrum_get_nwave(const cpl_table *table)
char GRAVI_DATAERR[50][10]
short gravi_sta_index(int gravi_input, cpl_table *optical_train_table, cpl_table *array_geometry_table)
Retrieve STA_INDEX corresponding to a given input.
int GRAVI_BASE_TEL[GRAVI_NBASE][2]
cpl_table * gravi_table_oi_create(int nwave, int nrow, const char *oi_name)
Create the oi table (oi_vis, oi_vis2, oi_t3)
int GRAVI_TRI_BASE[GRAVI_NBASE][2][2]