#include <complex.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include "mwa_hyperbeam.h"
#ifdef SINGLE
#define FLOAT float
#define CREAL crealf
#define CIMAG cimagf
#else
#define FLOAT double
#define CREAL creal
#define CIMAG cimag
#endif
void handle_hyperbeam_error(char file[], int line_num, const char function_name[]) {
int err_length = hb_last_error_length();
char *err = malloc(err_length * sizeof(char));
int err_status = hb_last_error_message(err, err_length);
if (err_status == -1) {
printf("Something really bad happened!\n");
exit(EXIT_FAILURE);
}
printf("File %s:%d: hyperbeam error in %s: %s\n", file, line_num, function_name, err);
exit(EXIT_FAILURE);
}
int main(int argc, char *argv[]) {
AnalyticBeam *beam;
char rts_style = 0; double *dipole_height_metres = NULL; uint8_t *bowties_per_row = NULL;
if (new_analytic_beam(rts_style, dipole_height_metres, bowties_per_row, &beam))
handle_hyperbeam_error(__FILE__, __LINE__, "new_analytic_beam");
unsigned delays[32] = {3, 2, 1, 0, 3, 2, 1, 0, 3, 2, 1, 0, 3, 2, 1, 0,
3, 2, 1, 0, 3, 2, 1, 0, 3, 2, 1, 0, 3, 2, 1, 0};
double dip_amps[32] = {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0, 1,
1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0};
int freqs_hz[3] = {150e6, 175e6, 200e6};
int num_freqs = 3;
int num_tiles = 2;
int num_amps = 16;
double latitude_rad = -0.4660608448386394;
int norm_to_zenith = 1;
AnalyticBeamGpu *gpu_beam;
if (new_gpu_analytic_beam(beam, delays, dip_amps, num_tiles, num_amps, &gpu_beam))
handle_hyperbeam_error(__FILE__, __LINE__, "new_gpu_analytic_beam");
int num_azzas = 1000000;
FLOAT *az = malloc(num_azzas * sizeof(FLOAT));
FLOAT *za = malloc(num_azzas * sizeof(FLOAT));
for (int i = 0; i < num_azzas; i++) {
az[i] = (-170.0 + i * 340.0 / num_azzas) * M_PI / 180.0;
za[i] = (10.0 + i * 70.0 / num_azzas) * M_PI / 180.0;
}
size_t num_unique_tiles = (size_t)get_num_unique_analytic_tiles(gpu_beam);
complex FLOAT *jones = malloc(num_unique_tiles * num_freqs * num_azzas * 8 * sizeof(FLOAT));
if (analytic_calc_jones_gpu(gpu_beam, num_azzas, az, za, num_freqs, freqs_hz, latitude_rad, norm_to_zenith,
(FLOAT *)jones))
handle_hyperbeam_error(__FILE__, __LINE__, "analytic_calc_jones_gpu");
printf("The first Jones matrix:\n");
printf("[[%+.8f%+.8fi,", CREAL(jones[0]), CIMAG(jones[0]));
printf(" %+.8f%+.8fi]\n", CREAL(jones[1]), CIMAG(jones[1]));
printf(" [%+.8f%+.8fi,", CREAL(jones[2]), CIMAG(jones[2]));
printf(" %+.8f%+.8fi]]\n", CREAL(jones[3]), CIMAG(jones[3]));
free(jones);
free_gpu_analytic_beam(gpu_beam);
free_analytic_beam(beam);
return EXIT_SUCCESS;
}