Skip to content
Snippets Groups Projects
Select Git revision
  • eda5263ee4a267769d0e14b3eedbc7503d12a0c4
  • master default protected
  • initialize_telescope_class_without_ms
  • rtd-ubuntu24-plucky
  • ast-1644-temp-fix
  • ast-1600-fix-beam-for-meerkat-ska-mid
  • readthedocs-c++17
  • mwa_python_wrapper
  • ast-1384-remove-element-index-argument
  • ast-1509-fix-polarization-orientation-in-gridded-response
  • add-test-for-squared-mueller-matrices
  • 1493-extend-python-bindings
  • ast-1493-implement-response-dishpoint-1
  • ast-1325-prototype-ska-beam-model-interface
  • ast-1416-oskar-ska-sdp-func-1
  • ast-1386-create-default-element
  • ast-1384-fix-sorted-frequencies-check-sphericalharmonicsresponse-1
  • ast-1111-add-vector-bindings
  • ast-973-add-test-for-lobes-coefficients
  • ast-645-add-beam-normalisation-mode-preapplied
  • disable-element-beam-1
  • v0.7.2
  • v0.7.1
  • v0.7.0
  • v0.6.2
  • v0.6.1
  • v0.6.0
  • v0.5.8
  • v0.5.7
  • v0.5.6
  • v0.5.5
  • v0.5.4
  • v0.5.3
  • v0.5.2
  • v0.5.1
  • v0.4.0
  • v0.3.1
  • v0.3.0
  • v0.2.0
  • v0.1.3
  • v0.1.2
41 results

tElementBeamCommon.h

Blame
  • Code owners
    Assign users and groups as approvers for specific file changes. Learn more.
    tElementBeamCommon.h 5.38 KiB
    #include <iostream>
    #include <complex>
    #include <vector>
    
    #include "beam-helper.h"
    
    void calculateElementBeams(everybeam::Station::Ptr& station,
                               std::vector<vector2r_t>& thetaPhiDirections,
                               size_t nr_antennas, unsigned int subgrid_size,
                               double frequency,
                               std::vector<std::complex<float>>& buffer) {
      typedef std::complex<float> Data[nr_antennas][subgrid_size][subgrid_size][4];
      Data* data_ptr = (Data*)buffer.data();
    
      auto elementResponse = station->get_element_response();
    
    #pragma omp parallel for
      for (size_t a = 0; a < nr_antennas; a++) {
        for (unsigned y = 0; y < subgrid_size; y++) {
          for (unsigned x = 0; x < subgrid_size; x++) {
            // Get theta, phi
            auto direction_thetaphi = thetaPhiDirections[y * subgrid_size + x];
            double theta = direction_thetaphi[0];
            double phi = direction_thetaphi[1];
    
            // Compute gain
            std::complex<double> gainMatrix[2][2] = {0.0};
            if (std::isfinite(theta) && std::isfinite(phi)) {
              elementResponse->Response(a, frequency, theta, phi, gainMatrix);
            }
    
            // Store gain
            std::complex<float>* antBufferPtr = (*data_ptr)[a][y][x];
            antBufferPtr[0] = gainMatrix[0][0];
            antBufferPtr[1] = gainMatrix[0][1];
            antBufferPtr[2] = gainMatrix[1][0];
            antBufferPtr[3] = gainMatrix[1][1];
          }
        }
      }
    }
    
    void calculateElementBeams(everybeam::Station::Ptr& station,
                               std::vector<vector3r_t>& itrfDirections,
                               size_t nr_antennas, unsigned int subgrid_size,
                               double time, double frequency,
                               std::vector<std::complex<float>>& buffer) {
      typedef std::complex<float> Data[nr_antennas][subgrid_size][subgrid_size][4];
      Data* data_ptr = (Data*)buffer.data();
    
    #pragma omp parallel for
      for (size_t a = 0; a < nr_antennas; a++) {
        for (unsigned y = 0; y < subgrid_size; y++) {
          for (unsigned x = 0; x < subgrid_size; x++) {
            // Get direction
            auto direction = itrfDirections[y * subgrid_size + x];
    
            // Compute gain
            matrix22c_t gainMatrix = {0.0};
            if (std::isfinite(direction[0])) {
              gainMatrix =
                  station->elementResponse(time, frequency, direction, a, true);
            }
    
            // Store gain
            std::complex<float>* antBufferPtr = (*data_ptr)[a][y][x];
            antBufferPtr[0] = gainMatrix[0][0];
            antBufferPtr[1] = gainMatrix[0][1];
            antBufferPtr[2] = gainMatrix[1][0];
            antBufferPtr[3] = gainMatrix[1][1];
          }
        }
      }
    }
    
    void run(everybeam::ElementResponseModel elementResponseModel, double frequency,
             std::string& input_filename, std::string& output_filename) {
      // Open measurement set
      std::cout << ">> Opening measurementset: " << input_filename << std::endl;
      casacore::MeasurementSet ms(input_filename);
    
      // Print frequency
      std::clog << "Frequency: " << frequency * 1e-6 << " Mhz" << std::endl;
    
      // Set number of stations to 1
      size_t nr_stations = 1;
    
      // Read number of timesteps
      size_t nr_timesteps = ms.nrow();
      std::clog << "Number of timesteps: " << nr_timesteps << std::endl;
    
      // Read observation time
      casacore::ScalarColumn<double> timeColumn(
          ms, ms.columnName(casacore::MSMainEnums::TIME));
      double currentTime = timeColumn(nr_timesteps / 2);
    
      // Read station
      size_t field_id = 0;
      size_t station_id = 0;
      auto station = readStation(ms, station_id, elementResponseModel);
      auto field_name = GetFieldName(ms, field_id);
      auto station_name = GetStationName(ms, station_id);
      auto nr_antennas = GetNrAntennas(ms, field_id);
      std::cout << "field: " << field_name << std::endl;
      std::cout << "station: " << station_name << std::endl;
      std::cout << "nr_antennas: " << nr_antennas << std::endl;
    
      // Compute RA and DEC of zenith at currentTime
      double zenithRA, zenithDec;
      GetRaDecZenith(station->position(), currentTime, zenithRA, zenithDec);
      std::clog << "RA:  " << zenithRA << std::endl;
      std::clog << "DEC: " << zenithDec << std::endl;
    
      // Imaging parameters
      float image_size = M_PI;   // in radians
      size_t subgrid_size = 32;  // in pixels
    
      // Compute element beams from theta, phi
      std::cout << ">>> Computing element beams theta, phi" << std::endl;
      std::vector<vector2r_t> thetaPhiDirections(subgrid_size * subgrid_size);
      GetThetaPhiDirectionsZenith(thetaPhiDirections.data(), subgrid_size);
      std::vector<std::complex<float>> beam_thetaphi(subgrid_size * subgrid_size *
                                                     4 * nr_antennas);
      calculateElementBeams(station, thetaPhiDirections, nr_antennas, subgrid_size,
                            frequency, beam_thetaphi);
    
    // Compute element beams from itrf coordinates
    // TODO: the Station::elementResponse method does not work properly
    #if 0
        std::cout << ">>> Computing element beams itrfs" << std::endl;
        std::vector<vector3r_t> itrfDirections(subgrid_size * subgrid_size);
        GetITRFDirections(itrfDirections.data(), subgrid_size, image_size, currentTime, zenithRA, zenithDec);
        std::vector<std::complex<float>> beam_itrf(subgrid_size*subgrid_size*4*nr_antennas);
        calculateElementBeams(station, itrfDirections, nr_antennas, subgrid_size, currentTime, frequency, beam_itrf);
    #endif
    
      // Store aterm
      StoreBeam(output_filename, beam_thetaphi.data(), nr_antennas, subgrid_size,
                subgrid_size);
    }