Open source Very Long Baseline Interferometry
OpenVLBI
Loading...
Searching...
No Matches
vlbi.h
1/* OpenVLBI - Open Source Very Long Baseline Interferometry
2* Copyright © 2017-2023 Ilia Platone
3*
4* This program is free software; you can redistribute it and/or modify
5* it under the terms of the GNU General Public License as published by
6* the Free Software Foundation; either version 2 of the License, or
7* (at your option) any later version.
8*
9* This program is distributed in the hope that it will be useful,
10* but WITHOUT ANY WARRANTY; without even the implied warranty of
11* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12* GNU General Public License for more details.
13*
14* You should have received a copy of the GNU General Public License along
15* with this program; if not, write to the Free Software Foundation, Inc.,
16* 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
17*/
18
19#ifndef _VLBI_H
20#define _VLBI_H
21
22#ifdef __cplusplus
23extern "C" {
24#endif
25
26#ifndef DLL_EXPORT
27#ifdef _WIN32
28#define DLL_EXPORT __declspec(dllexport)
29#else
30#define DLL_EXPORT extern
31#endif
32#endif
33
34#include <stdio.h>
35#include <stdlib.h>
36#include <math.h>
37#include <time.h>
38#include <errno.h>
39#include <unistd.h>
40#include <string.h>
41#include <math.h>
42#include <float.h>
43#include <libgen.h>
44#include <sys/types.h>
45#include <assert.h>
46#include <pthread.h>
47#include <dsp.h>
48
283typedef struct {
287 double *Location;
289 int Geo;
293 char *Name;
295 int Index;
296} vlbi_node;
297
299typedef struct {
305 double *Target;
307 double Ra;
309 double Dec;
311 double *baseline;
313 double u;
315 double v;
317 double delay;
327 char *Name;
331
342typedef double(* vlbi_func2_t)(double, double);
343
345typedef void* vlbi_context;
346
348typedef struct timespec timespec_t;
355#define VLBI_BASELINE_NAME_SIZE DSP_NAME_SIZE*2+1
356
357#ifndef Min
359#define Min(a,b) \
360 ({ __typeof (a) _a = (a); \
361 __typeof (b) _b = (b); \
362 _a < _b ? _a : _b; })
363#endif
364#ifndef Max
366#define Max(a,b) \
367 ({ __typeof (a) _a = (a); \
368 __typeof (b) _b = (b); \
369 _a > _b ? _a : _b; })
370#endif
371#ifndef Log
373#define Log(a,b) \
374( log(a) / log(b) )
375#endif
376
377#ifndef hz2rad
379#define hz2rad(hz) (2.0*PI*hz)
380#endif
381
382#ifndef sinpsin
384#define sinpsin(s1, s2) (2.0*sin((asin(s1)+asin(s2))/2.0)*cos((asin(s1)-asin(s2))/2.0))
385#endif
386
387#ifndef sinmsin
389#define sinmsin(s1, s2) (2.0*cos((asin(s1)+asin(s2))/2.0)*sin((asin(s1)-asin(s2))/2.0))
390#endif
391
392#ifndef cospcos
394#define cospcos(c1, c2) (2.0*cos((acos(s1)+acos(s2))/2.0)*cos((acos(s1)-acos(s2))/2.0))
395#endif
396
397#ifndef cosmcos
399#define cosmcos(c1, c2) (2.0*sin((acos(s1)+acos(s2))/2.0)*sin((acos(s1)-acos(s2))/2.0))
400#endif
401
402#ifndef sinxsin
404#define sinxsin(s1, s2) ((cos(asin(s1)-asin(s2))-cos(asin(s1)+asin(s2)))/2.0)
405#endif
406
407#ifndef cosxcos
409#define cosxcos(c1, c2) ((cos(acos(s1)+acos(s2))+cos(acos(s1)-acos(s2)))/2.0)
410#endif
411
412#ifndef sinxcos
414#define sinxcos(s, c) ((sin(asin(s)+acos(c))+sin(asin(s)-acos(c)))/2.0)
415#endif
416
417#ifndef sin2cos
419#define sin2cos(s) cos(asin(s))
420#endif
421
422#ifndef cos2sin
424#define cos2sin(c) sin(acos(c))
425#endif
426
427#ifndef VLBI_VERSION_STRING
429#define VLBI_VERSION_STRING "3.0.2"
430#endif
431
432#ifndef VLBI_CATALOG_PATH
434#define VLBI_CATALOG_PATH "/usr/share/OpenVLBI/cat/index.txt"
435#endif
436
437#ifndef CIRCLE_DEG
439#define CIRCLE_DEG 360.0
440#endif
441#ifndef CIRCLE_AM
443#define CIRCLE_AM (CIRCLE_DEG * 60.0)
444#endif
445#ifndef CIRCLE_AS
447#define CIRCLE_AS (CIRCLE_AM * 60.0)
448#endif
449#ifndef RAD_AS
451#define RAD_AS (CIRCLE_AS/(PI*2.0))
452#endif
453#ifndef ONE_SECOND_TICKS
455#define ONE_SECOND_TICKS 100000000
456#endif
457#ifndef ONE_MILLISECOND_TICKS
459#define ONE_MILLISECOND_TICKS 100000
460#endif
461#ifndef ONE_MICROSECOND_TICKS
463#define ONE_MICROSECOND_TICKS 100
464#endif
465#ifndef SOLAR_DAY
467#define SOLAR_DAY 86400
468#endif
469#ifndef SIDEREAL_DAY
471#define SIDEREAL_DAY 86164.0905
472#endif
473#ifndef TRACKRATE_SIDEREAL
475#define TRACKRATE_SIDEREAL (CIRCLE_AS / SIDEREAL_DAY)
476#endif
477#ifndef TRACKRATE_SOLAR
479#define TRACKRATE_SOLAR (CIRCLE_AS / SOLAR_DAY)
480#endif
481#ifndef TRACKRATE_LUNAR
483#define TRACKRATE_LUNAR 14.511415
484#endif
485#ifndef EARTHRADIUSEQUATORIAL
487#define EARTHRADIUSEQUATORIAL 6378137.0
488#endif
489#ifndef EARTHRADIUSPOLAR
491#define EARTHRADIUSPOLAR 6356752.0
492#endif
493#ifndef EARTHRADIUSMEAN
495#define EARTHRADIUSMEAN 6372797.0
496#endif
497#ifndef AVOGADRO
499#define AVOGADRO 6.02214076E+23
500#endif
501#ifndef EULER
503#define EULER 2.71828182845904523536028747135266249775724709369995
504#endif
505#ifndef PLANK
507#define PLANK 6.62607015E-34
508#endif
509#ifndef BOLTSMANN
511#define BOLTSMANN 1.380649E-23
512#endif
513#ifndef STEPHAN_BOLTSMANN
515#define STEPHAN_BOLTSMANN (2.0*pow(PI, 5)*pow(BOLTSMANN, 2)/(pow(LIGHTSPEED, 2)*15*pow(PLANK, 3)))
516#endif
517#ifndef GAS_R
519#define GAS_R (BOLTSMANN * AVOGADRO)
520#endif
521#ifndef ROOT2
523#define ROOT2 1.41421356237309504880168872420969807856967187537694
524#endif
525#ifndef PI
527#define PI 3.14159265358979323846
528#endif
529#ifndef AIRY
531#define AIRY 1.21966
532#endif
533#ifndef LIGHTSPEED
535#define LIGHTSPEED 299792458.0
536#endif
537#ifndef J2000
539#define J2000 2451545.0
540#endif
541#ifndef GAMMAJ2000
543#define GAMMAJ2000 18.6971378528
544#endif
545#ifndef ELECTRON
547#define ELECTRON 1.602176634E-19
548#endif
549#ifndef CANDLE
551#define CANDLE 0.683
552#endif
553#ifndef ASTRONOMICALUNIT
555#define ASTRONOMICALUNIT 1.495978707E+11
556#endif
557#ifndef PARSEC
559#define PARSEC (ASTRONOMICALUNIT/sin(PI*2.0/CIRCLE_AS))
560#endif
561#ifndef LY
563#define LY (LIGHTSPEED * SIDEREAL_DAY * 365.0)
564#endif
565#ifndef AU2M
571#define AU2M(au) (au * ASTRONOMICALUNIT)
572#endif
573#ifndef PARSEC2M
579#define PARSEC2M(parsec) (parsec * PARSEC)
580#endif
581#ifndef LY2M
587#define LY2M(ly) (ly * LY)
588#endif
589#ifndef M2AU
595#define M2AU(m) (m / ASTRONOMICALUNIT)
596#endif
597#ifndef M2PARSEC
603#define M2PARSEC(m) (m / PARSEC)
604#endif
605#ifndef M2LY
611#define M2LY(m) (m / LY)
612#endif
613#ifndef RAD2AS
619#define RAD2AS(rad) (rad * RAD_AS)
620#endif
621#ifndef AS2RAD
627#define AS2RAD(as) (as / RAD_AS);
628#endif
629
637inline double vlbi_default_delegate(double x, double y) {
638 return x*y;
639}
640
649inline double vlbi_magnitude_delegate(double x, double y) {
650 return sqrt(pow(x, 2)+pow(y, 2));
651}
652
661inline double vlbi_phase_delegate(double x, double y) {
662 double mag = sqrt(pow(x, 2)+pow(y, 2));
663 double rad = 0.0;
664 if(mag > 0.0) {
665 rad = acos (y / (mag > 0.0 ? mag : 1.0));
666 if(x < 0 && rad != 0)
667 rad = PI*2-rad;
668 }
669 return rad;
670}
671
679inline double vlbi_magnitude_correlator_delegate(double x, double y) {
680 return x*y;
681}
682
690inline double vlbi_phase_correlator_delegate(double x, double y) {
691 return sinxsin(x, y);
692}
693
705DLL_EXPORT unsigned long int vlbi_max_threads(unsigned long value);
706
711DLL_EXPORT const char *vlbi_get_version(void);
712
717DLL_EXPORT vlbi_context vlbi_init(void);
718
723DLL_EXPORT void vlbi_exit(vlbi_context ctx);
724
738DLL_EXPORT void vlbi_add_node(vlbi_context ctx, dsp_stream_p Stream, const char *name, int geographic_coordinates);
739
746DLL_EXPORT void vlbi_copy_node(void *ctx, const char *name, const char *node);
747
754DLL_EXPORT dsp_stream_p vlbi_get_node(void *ctx, const char *name);
755
762DLL_EXPORT int vlbi_has_node(void *ctx, const char *name);
763
769DLL_EXPORT void vlbi_del_node(vlbi_context ctx, const char *name);
770
777DLL_EXPORT int vlbi_get_nodes(void *ctx, vlbi_node** nodes);
778
786DLL_EXPORT void vlbi_add_node_from_fits(void *ctx, char *filename, const char *name, int geo);
787
795DLL_EXPORT void vlbi_add_nodes_from_sdfits(void *ctx, char *filename, const char *name, int geo);
796
804DLL_EXPORT void vlbi_filter_lp_node(void *ctx, const char *name, const char *node, double radians);
805
813DLL_EXPORT void vlbi_filter_hp_node(void *ctx, const char *name, const char *node, double radians);
814
823DLL_EXPORT void vlbi_filter_bp_node(void *ctx, const char *name, const char *node, double lo_radians, double hi_radians);
824
833DLL_EXPORT void vlbi_filter_br_node(void *ctx, const char *name, const char *node, double lo_radians, double hi_radians);
834
846DLL_EXPORT void vlbi_set_correlation_order(void *ctx, int order);
847
854DLL_EXPORT int vlbi_get_baselines(void *ctx, vlbi_baseline** baselines);
855
861DLL_EXPORT void vlbi_free_baselines(vlbi_baseline** baselines, int num_baselines);
862
873DLL_EXPORT void vlbi_set_baseline_buffer(void *ctx, const char **nodes, complex_t *buffer, int len);
874
881DLL_EXPORT dsp_stream_p vlbi_get_baseline_stream(void *ctx, const char**nodes);
882
888DLL_EXPORT void vlbi_unlock_baseline(void *ctx, const char**nodes);
889
899DLL_EXPORT void vlbi_set_baseline_stream(void *ctx, const char**nodes, dsp_stream_p stream);
900
908DLL_EXPORT void vlbi_set_location(void *ctx, double lat, double lon, double el);
909
920DLL_EXPORT double vlbi_get_offset(vlbi_context ctx, double J2000Time, const char* node, double Ra, double Dec, double Distance);
921
942DLL_EXPORT void vlbi_get_uv_plot(void *ctx, const char *name, int u, int v, double *target, double freq, double sr, int nodelay, int moving_baseline, vlbi_func2_t delegate, int *interrupt);
943
950DLL_EXPORT void vlbi_add_model(vlbi_context ctx, dsp_stream_p Stream, const char *name);
951
958DLL_EXPORT void vlbi_copy_model(void *ctx, const char *name, const char *node);
959
965DLL_EXPORT void vlbi_del_model(vlbi_context ctx, const char *name);
966
973DLL_EXPORT int vlbi_get_models(void *ctx, dsp_stream_p** models);
974
981DLL_EXPORT dsp_stream_p vlbi_get_model(void *ctx, const char *name);
982
989DLL_EXPORT int vlbi_has_model(void *ctx, const char *name);
990
998DLL_EXPORT void vlbi_get_ifft(vlbi_context ctx, const char *name, const char *magnitude, const char *phase);
999
1007DLL_EXPORT void vlbi_get_fft(vlbi_context ctx, const char *model, const char *magnitude, const char *phase);
1008
1016DLL_EXPORT void vlbi_apply_mask(vlbi_context ctx, const char *name, const char *model, const char *mask);
1017
1025DLL_EXPORT void vlbi_apply_convolution_matrix(vlbi_context ctx, const char *name, const char *model, const char *matrix);
1026
1034DLL_EXPORT void vlbi_stack_models(vlbi_context ctx, const char *name, const char *model1, const char *model2);
1035
1043DLL_EXPORT void vlbi_diff_models(vlbi_context ctx, const char *name, const char *model1, const char *model2);
1044
1050DLL_EXPORT void vlbi_shift(vlbi_context ctx, const char *name);
1051
1058DLL_EXPORT void vlbi_add_model_from_png(void *ctx, char *filename, const char *name);
1059
1066DLL_EXPORT void vlbi_add_model_from_jpeg(void *ctx, char *filename, const char *name);
1067
1074DLL_EXPORT void vlbi_add_model_from_fits(void *ctx, char *filename, const char *name);
1075
1082DLL_EXPORT void vlbi_get_model_to_png(void *ctx, char *filename, const char *name);
1083
1090DLL_EXPORT void vlbi_get_model_to_jpeg(void *ctx, char *filename, const char *name);
1091
1098DLL_EXPORT void vlbi_get_model_to_fits(void *ctx, char *filename, const char *name);
1099
1113DLL_EXPORT double* vlbi_matrix_calc_baseline_center(double *loc1, double *loc2);
1114
1122DLL_EXPORT double* vlbi_matrix_calc_3d_projection(double alt, double az, double *baseline);
1123
1130DLL_EXPORT double* vlbi_matrix_calc_parametric_projection(double *target, double *baseline);
1131
1138DLL_EXPORT double* vlbi_matrix_calc_uv_coordinates(double *proj, double wavelength);
1139
1145DLL_EXPORT double* vlbi_matrix_calc_location(double *loc);
1146
1153DLL_EXPORT double vlbi_matrix_estimate_geocentric_elevation(double latitude, double sea_level_elevation);
1154
1160DLL_EXPORT double vlbi_matrix_estimate_resolution_zero(double frequency);
1161
1168DLL_EXPORT double vlbi_matrix_estimate_resolution(double resolution0, double baseline);
1169
1187DLL_EXPORT timespec_t vlbi_time_mktimespec(int year, int month, int dom, int hour, int minute, int second, long nanosecond);
1188
1195
1201DLL_EXPORT double vlbi_time_J2000time_to_lst(double secs_since_J2000, double Long);
1202
1208DLL_EXPORT timespec_t vlbi_time_string_to_timespec(const char *time);
1209
1215DLL_EXPORT timespec_t vlbi_time_J2000time_to_timespec(double secs_since_J2000);
1216
1228DLL_EXPORT double vlbi_astro_mean_speed(double speed);
1229
1240DLL_EXPORT void vlbi_astro_alt_az_from_ra_dec(double J2000time, double Ra, double Dec, double Lat, double Long, double* Alt, double *Az);
1241
1248DLL_EXPORT double vlbi_astro_get_local_hour_angle(double local_sideral_time, double ra);
1249
1258DLL_EXPORT void vlbi_astro_get_alt_az_coordinates(double hour_angle, double dec, double latitude, double* alt, double *az);
1259
1265DLL_EXPORT dsp_stream_p vlbi_astro_load_spectrum(char *filename);
1266
1274DLL_EXPORT int vlbi_astro_load_spectra_catalog(char *path, dsp_stream_p **catalog, int *catalog_size);
1275
1283
1289DLL_EXPORT void vlbi_astro_save_spectrum(dsp_stream_p stream, char *filename);
1290
1296DLL_EXPORT void vlbi_astro_scan_spectrum(dsp_stream_p stream, int sample_size);
1297
1307DLL_EXPORT dsp_align_info vlbi_astro_align_spectra(dsp_stream_p spectrum, dsp_stream_p catalog, int max_lines, double decimals, double min_score);
1308
1316DLL_EXPORT double vlbi_astro_diff_spectra(dsp_stream_p spectrum0, dsp_stream_p spectrum, double decades);
1317
1325DLL_EXPORT double vlbi_astro_flux_ratio(double flux0, double flux, double delta_spectrum);
1326
1333DLL_EXPORT double vlbi_astro_estimate_brightness_temperature(double wavelength, double flux);
1334
1341DLL_EXPORT double vlbi_astro_estimate_temperature(double wavelength, double flux);
1342
1349DLL_EXPORT double vlbi_astro_estimate_flux(double wavelength, double temperature);
1350
1357DLL_EXPORT double vlbi_astro_estimate_temperature_ratio(double rad_ratio, double flux_ratio);
1358
1365DLL_EXPORT double vlbi_astro_estimate_size_ratio(double luminosity_ratio, double temperature_ratio);
1366
1373DLL_EXPORT double vlbi_astro_estimate_luminosity_ratio(double size_ratio, double flux_ratio);
1374
1381DLL_EXPORT double vlbi_astro_estimate_distance_ratio(double luminosity_ratio, double flux_ratio);
1382
1389DLL_EXPORT double vlbi_astro_estimate_distance_parallax(double rad, double baseline);
1390
1397DLL_EXPORT double vlbi_astro_estimate_redshift(double wavelength0, int wavelength);
1398
1405DLL_EXPORT double vlbi_astro_estimate_size_transient(double transient_object_velocity, double transit_time);
1406
1413DLL_EXPORT double vlbi_astro_redshift_adjust(double distance, double redshift);
1414
1420DLL_EXPORT dsp_stream_p vlbi_file_read_fits(char *filename);
1421
1428DLL_EXPORT dsp_stream_p *vlbi_file_read_sdfits(char * filename, long *n);
1429
1435#ifdef __cplusplus
1436}
1437#endif
1438
1439#endif //_VLBI_H
1440
double complex_t[2]
complex type
Definition dsp.h:77
DLL_EXPORT double vlbi_astro_estimate_temperature(double wavelength, double flux)
Estimate the physical temperature with a given flux.
DLL_EXPORT dsp_stream_p vlbi_astro_create_reference_catalog(dsp_stream_p *catalog, int catalog_size)
Create a dsp_stream_p containing all the spectral lines of all elements of a catalog.
DLL_EXPORT double vlbi_astro_diff_spectra(dsp_stream_p spectrum0, dsp_stream_p spectrum, double decades)
Compare a spectrum to a reference one.
DLL_EXPORT dsp_stream_p * vlbi_file_read_sdfits(char *filename, long *n)
Returns the distance of a far object adjusted with its measured redshift.
DLL_EXPORT double vlbi_astro_estimate_distance_ratio(double luminosity_ratio, double flux_ratio)
Returns the distance ratio of two celestial object.
DLL_EXPORT dsp_align_info vlbi_astro_align_spectra(dsp_stream_p spectrum, dsp_stream_p catalog, int max_lines, double decimals, double min_score)
Align a spectrum to a reference catalog.
DLL_EXPORT double vlbi_astro_estimate_flux(double wavelength, double temperature)
Estimate the flux from a physical temperature.
DLL_EXPORT double vlbi_astro_estimate_distance_parallax(double rad, double baseline)
Returns the distance of an object from its parallax.
DLL_EXPORT dsp_stream_p vlbi_astro_load_spectrum(char *filename)
Load a spectrum file.
DLL_EXPORT int vlbi_astro_load_spectra_catalog(char *path, dsp_stream_p **catalog, int *catalog_size)
Load a spectrum file catalog.
DLL_EXPORT double vlbi_astro_estimate_brightness_temperature(double wavelength, double flux)
Estimate the brightness temperature from a flux.
DLL_EXPORT double vlbi_astro_flux_ratio(double flux0, double flux, double delta_spectrum)
Returns the flux ratio of two objects.
DLL_EXPORT void vlbi_astro_alt_az_from_ra_dec(double J2000time, double Ra, double Dec, double Lat, double Long, double *Alt, double *Az)
Obtain the altitude and azimuth coordinate of a celestial coordinate at a specific time.
DLL_EXPORT void vlbi_astro_scan_spectrum(dsp_stream_p stream, int sample_size)
Scan a dsp_stream_p and detect the relevant spectral lines.
DLL_EXPORT void vlbi_astro_save_spectrum(dsp_stream_p stream, char *filename)
Save a spectrum stream into a file.
DLL_EXPORT void vlbi_astro_get_alt_az_coordinates(double hour_angle, double dec, double latitude, double *alt, double *az)
Returns alt-azimuth coordinates of an object.
DLL_EXPORT double vlbi_astro_estimate_luminosity_ratio(double size_ratio, double flux_ratio)
Returns the luminosity ratio between two celestial objects.
DLL_EXPORT double vlbi_astro_redshift_adjust(double distance, double redshift)
Returns the distance of a far object adjusted with its measured redshift.
DLL_EXPORT double vlbi_astro_estimate_temperature_ratio(double rad_ratio, double flux_ratio)
Returns the temperature ratio of two objects.
DLL_EXPORT double vlbi_astro_estimate_redshift(double wavelength0, int wavelength)
Returns the redshift of an object.
DLL_EXPORT dsp_stream_p vlbi_file_read_fits(char *filename)
Returns the distance of a far object adjusted with its measured redshift.
DLL_EXPORT double vlbi_astro_estimate_size_transient(double transient_object_velocity, double transit_time)
Returns the size of an object by the observation of a transient body.
DLL_EXPORT double vlbi_astro_estimate_size_ratio(double luminosity_ratio, double temperature_ratio)
Returns the size ratio of two objects.
DLL_EXPORT double vlbi_astro_get_local_hour_angle(double local_sideral_time, double ra)
Returns local hour angle of an object.
DLL_EXPORT double vlbi_astro_mean_speed(double speed)
Obtain or set the reference constant speed of the radiation to measure.
DLL_EXPORT double vlbi_get_offset(vlbi_context ctx, double J2000Time, const char *node, double Ra, double Dec, double Distance)
Get the offset of a single node to the farest node to the target.
DLL_EXPORT void vlbi_unlock_baseline(void *ctx, const char **nodes)
Unlock the baseline and get visibility from its nodes correlations.
DLL_EXPORT dsp_stream_p vlbi_get_baseline_stream(void *ctx, const char **nodes)
Obtain the baseline dsp_stream structure containing the complex visibility data.
DLL_EXPORT void vlbi_set_baseline_stream(void *ctx, const char **nodes, dsp_stream_p stream)
Set the baseline dsp_stream structure containing the complex visibility data. This function locks thi...
DLL_EXPORT void vlbi_set_correlation_order(void *ctx, int order)
Set the correlation order the current OpenVLBI context.
DLL_EXPORT void vlbi_set_baseline_buffer(void *ctx, const char **nodes, complex_t *buffer, int len)
Fill the buffer of a single baseline with complex visibility data. This function locks this baeline a...
DLL_EXPORT int vlbi_get_baselines(void *ctx, vlbi_baseline **baselines)
List all baselines of the current OpenVLBI context.
DLL_EXPORT void vlbi_set_location(void *ctx, double lat, double lon, double el)
Set the location of the reference station.
DLL_EXPORT void vlbi_free_baselines(vlbi_baseline **baselines, int num_baselines)
Free a vlbi_baseline object previously allocated by vlbi_get_baselines.
#define sinxsin(s1, s2)
Multiply a sine to a sine.
Definition vlbi.h:404
double vlbi_default_delegate(double x, double y)
A placeholder delegate that simply multiplies the values received from vlbi_get_uv_plot.
Definition vlbi.h:637
#define PI
Our PI constant.
Definition vlbi.h:527
double vlbi_magnitude_correlator_delegate(double x, double y)
A magnitude correlator delegate for vlbi_get_uv_plot.
Definition vlbi.h:679
double vlbi_phase_correlator_delegate(double x, double y)
A phase correlator delegate for vlbi_get_uv_plot.
Definition vlbi.h:690
double vlbi_phase_delegate(double x, double y)
A phase calculator delegate for vlbi_get_uv_plot.
Definition vlbi.h:661
double vlbi_magnitude_delegate(double x, double y)
A magnitude calculator delegate for vlbi_get_uv_plot.
Definition vlbi.h:649
DLL_EXPORT const char * vlbi_get_version(void)
Print the current version of OpenVLBI.
DLL_EXPORT void vlbi_exit(vlbi_context ctx)
Close a OpenVLBI instance.
DLL_EXPORT unsigned long int vlbi_max_threads(unsigned long value)
get/set the maximum number of threads allowed
DLL_EXPORT vlbi_context vlbi_init(void)
Initialize a OpenVLBI instance.
DLL_EXPORT double * vlbi_matrix_calc_baseline_center(double *loc1, double *loc2)
Return The baseline center in geographic coordinates.
DLL_EXPORT double * vlbi_matrix_calc_uv_coordinates(double *proj, double wavelength)
Return The UV coordinates of the current observation.
DLL_EXPORT double * vlbi_matrix_calc_3d_projection(double alt, double az, double *baseline)
Return The 3d projection of the current observation.
DLL_EXPORT double * vlbi_matrix_calc_parametric_projection(double *target, double *baseline)
Return The parametric projection of the current observation.
DLL_EXPORT double * vlbi_matrix_calc_location(double *loc)
Convert geographic location into xyz location.
DLL_EXPORT double vlbi_matrix_estimate_geocentric_elevation(double latitude, double sea_level_elevation)
Returns an estimation of the actual geocentric elevation.
DLL_EXPORT double vlbi_matrix_estimate_resolution_zero(double frequency)
Estimate the angular resolution of a 1 meter baseline at a given frequency.
DLL_EXPORT double vlbi_matrix_estimate_resolution(double resolution0, double baseline)
Estimate Signal to noise ratio after a given integration time.
DLL_EXPORT void vlbi_add_model_from_fits(void *ctx, char *filename, const char *name)
Add a model from a fits file.
DLL_EXPORT void vlbi_get_model_to_fits(void *ctx, char *filename, const char *name)
Write a model to a fits file.
DLL_EXPORT void vlbi_add_model(vlbi_context ctx, dsp_stream_p Stream, const char *name)
Add a model into the current OpenVLBI context.
DLL_EXPORT int vlbi_get_models(void *ctx, dsp_stream_p **models)
List all models of the current OpenVLBI context.
DLL_EXPORT void vlbi_get_model_to_jpeg(void *ctx, char *filename, const char *name)
Write a model to a jpeg file.
DLL_EXPORT void vlbi_stack_models(vlbi_context ctx, const char *name, const char *model1, const char *model2)
Stack two models into a new one.
DLL_EXPORT void vlbi_diff_models(vlbi_context ctx, const char *name, const char *model1, const char *model2)
Diff two models into a new one.
DLL_EXPORT int vlbi_has_model(void *ctx, const char *name)
Determine if a model is present into the current OpenVLBI context.
DLL_EXPORT void vlbi_add_model_from_png(void *ctx, char *filename, const char *name)
Add a model from a png file.
DLL_EXPORT void vlbi_copy_model(void *ctx, const char *name, const char *node)
Copy a model into a new one.
DLL_EXPORT void vlbi_shift(vlbi_context ctx, const char *name)
Shift a model by its dimensions.
DLL_EXPORT void vlbi_apply_convolution_matrix(vlbi_context ctx, const char *name, const char *model, const char *matrix)
Convolute a model with a convolution matrix.
DLL_EXPORT void vlbi_get_ifft(vlbi_context ctx, const char *name, const char *magnitude, const char *phase)
Save into name an inverse fourier transform of the uv plot using its current magnitude and phase comp...
DLL_EXPORT void vlbi_get_fft(vlbi_context ctx, const char *model, const char *magnitude, const char *phase)
Get the fourier transform of the given model and save its magnitude and phase components into two new...
DLL_EXPORT dsp_stream_p vlbi_get_model(void *ctx, const char *name)
Get a single model from the current OpenVLBI context.
DLL_EXPORT void vlbi_del_model(vlbi_context ctx, const char *name)
Remove a model from the current OpenVLBI context.
DLL_EXPORT void vlbi_add_model_from_jpeg(void *ctx, char *filename, const char *name)
Add a model from a jpeg file.
DLL_EXPORT void vlbi_apply_mask(vlbi_context ctx, const char *name, const char *model, const char *mask)
Mask the stream with the content of the mask stream, by multiplication of each element.
DLL_EXPORT void vlbi_get_model_to_png(void *ctx, char *filename, const char *name)
Write a model to a png file.
DLL_EXPORT void vlbi_get_uv_plot(void *ctx, const char *name, int u, int v, double *target, double freq, double sr, int nodelay, int moving_baseline, vlbi_func2_t delegate, int *interrupt)
Fill a fourier plane with an aperture synthesis projection of the baselines during the integration ti...
DLL_EXPORT dsp_stream_p vlbi_get_node(void *ctx, const char *name)
Get a stream from the current OpenVLBI context.
DLL_EXPORT void vlbi_copy_node(void *ctx, const char *name, const char *node)
Copy a node into a new one.
DLL_EXPORT void vlbi_add_node(vlbi_context ctx, dsp_stream_p Stream, const char *name, int geographic_coordinates)
Add a stream into the current OpenVLBI context.
DLL_EXPORT void vlbi_add_node_from_fits(void *ctx, char *filename, const char *name, int geo)
Add a node from a 2d image fits file.
DLL_EXPORT void vlbi_add_nodes_from_sdfits(void *ctx, char *filename, const char *name, int geo)
Add nodes from each row of a single dish fits -SDFITS- file.
DLL_EXPORT int vlbi_has_node(void *ctx, const char *name)
Determine if a node is present into the current OpenVLBI context.
DLL_EXPORT void vlbi_filter_bp_node(void *ctx, const char *name, const char *node, double lo_radians, double hi_radians)
Apply a band pass filter on the node buffer.
DLL_EXPORT void vlbi_del_node(vlbi_context ctx, const char *name)
Remove a stream from the current OpenVLBI context.
DLL_EXPORT void vlbi_filter_hp_node(void *ctx, const char *name, const char *node, double radians)
Apply a high pass filter on the node buffer.
DLL_EXPORT int vlbi_get_nodes(void *ctx, vlbi_node **nodes)
List all nodes of the current OpenVLBI context.
DLL_EXPORT void vlbi_filter_lp_node(void *ctx, const char *name, const char *node, double radians)
Apply a low pass filter on the node buffer.
DLL_EXPORT void vlbi_filter_br_node(void *ctx, const char *name, const char *node, double lo_radians, double hi_radians)
Apply a band reject filter on the node buffer.
DLL_EXPORT timespec_t vlbi_time_string_to_timespec(const char *time)
Obtain a timespec struct containing the date and time specified by a time string.
DLL_EXPORT timespec_t vlbi_time_J2000time_to_timespec(double secs_since_J2000)
Obtain a timespec struct containing the date and time specified by a J2000 time.
DLL_EXPORT timespec_t vlbi_time_mktimespec(int year, int month, int dom, int hour, int minute, int second, long nanosecond)
Obtain a timespec struct containing the date and time specified.
DLL_EXPORT double vlbi_time_timespec_to_J2000time(timespec_t tp)
Convert a timespec into J2000 time.
DLL_EXPORT double vlbi_time_J2000time_to_lst(double secs_since_J2000, double Long)
Obtain the local sidereal time at an arbitrary moment and location.
double(* vlbi_func2_t)(double, double)
The delegate function type to pass to vlbi_plot_uv_plane.
Definition vlbi.h:342
void * vlbi_context
the OpenVLBI context object type
Definition vlbi.h:345
struct timespec timespec_t
Definition of the timespec_t in a C type, just for convenience.
Definition vlbi.h:348
Alignment informations needed.
Definition dsp.h:294
Contains a set of informations and data relative to a buffer and how to use it.
Definition dsp.h:387
A single baseline, returned by vlbi_get_baselines into an array.
Definition vlbi.h:299
int relative
Whether this baseline's nodes use coordinates relative to the context reference station.
Definition vlbi.h:301
char * Name
Baseline's name.
Definition vlbi.h:327
double delay
Current delay in seconds.
Definition vlbi.h:317
double * baseline
The baseline's 3d sizes in meters.
Definition vlbi.h:311
double u
Current u perspective location in pixels.
Definition vlbi.h:313
int locked
Whether this baseline's nodes use coordinates relative to the context reference station.
Definition vlbi.h:303
double WaveLength
Wavelength observed in meters.
Definition vlbi.h:319
dsp_stream_p Stream
Baseline's DSP stream.
Definition vlbi.h:329
double * Target
The baseline's current celestial target array.
Definition vlbi.h:305
vlbi_node * Nodes
Nodes.
Definition vlbi.h:325
int nodes_count
Number of nodes.
Definition vlbi.h:323
double Ra
The baseline's current celestial target right ascension.
Definition vlbi.h:307
double SampleRate
Samples per second.
Definition vlbi.h:321
double Dec
The baseline's current celestial target declination.
Definition vlbi.h:309
double v
Current v perspective location in pixels.
Definition vlbi.h:315
A single node, returned by vlbi_get_nodes into an array.
Definition vlbi.h:283
int Geo
Whether this node uses geographic coordinates.
Definition vlbi.h:289
double * GeographicLocation
Geographic coordinates array (latitude, longitude, elevation)
Definition vlbi.h:285
dsp_stream_p Stream
Node's DSP stream.
Definition vlbi.h:291
int Index
Node's index - oldest zero.
Definition vlbi.h:295
double * Location
Current location array in meters.
Definition vlbi.h:287
char * Name
Node's name.
Definition vlbi.h:293