Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Implemented new definition of (u,v) points #171

Merged
merged 1 commit into from
Sep 21, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 5 additions & 8 deletions stix/idl/processing/imaging/stx_em.pro
Original file line number Diff line number Diff line change
Expand Up @@ -38,6 +38,7 @@
; June 2022, Massa P., 'aux_data' added
; August 2022, Massa P., made it compatible with the up-to-date imaging software and added backrgound correction
; in the EM iterative scheme
; July 2023, Massa P., made it compatible with the new definition of (u,v)-points (see stx_uv_points)
;
;CONTACT: massa.p@dima.unige.it

Expand Down Expand Up @@ -104,15 +105,11 @@ end
phase_corr = grid_phase_corr + ad_hoc_phase_corr + phase_factor
phase_corr *= !dtor

;;**************** Giordano's (u,v) points
;;**************** Define (u,v) points

subc_str = stx_construct_subcollimator()

uv = stx_uv_points_giordano()
u = -uv.u * subc_str.phase
v = -uv.v * subc_str.phase
u = u[subc_index]
v = v[subc_index]
uv_data = stx_uv_points(subc_index)
u = uv_data.u
v = uv_data.v

;;**************** Transmission matrix

Expand Down
13 changes: 2 additions & 11 deletions stix/idl/processing/imaging/stx_estimate_flare_location.pro
Original file line number Diff line number Diff line change
Expand Up @@ -38,6 +38,7 @@
; (default, 0)
;
; HISTORY: October 2022, Massa P. (WKU), initial release
; July 2023, Massa P., made it compatible with the new definition of (u,v)-points (see stx_uv_points)
;
; CONTACT:
; paolo.massa@wku.edu
Expand All @@ -62,20 +63,10 @@ pixel = rsun * 2.6 / imsize ; 2.6 is chosen arbitrarily so that field of view Ba
;;******* Construct the visibility structure
mapcenter_stix = stx_hpc2stx_coord(mapcenter, aux_data)

vis = stx_construct_visibility(path_sci_file, time_range, energy_range, mapcenter_stix, $
vis = stx_construct_calibrated_visibility(path_sci_file, time_range, energy_range, mapcenter_stix, $
path_bkg_file=path_bkg_file,subc_index=subc_index, /silent, $
_extra=extra)

;;******* Use Giordano's (u,v) points: no need to perform projection correction (see Giordano et al., 2015)
subc_str = stx_construct_subcollimator()
uv = stx_uv_points_giordano()
u = -uv.u * subc_str.phase
v = -uv.v * subc_str.phase
vis.u = u[subc_index]
vis.v = v[subc_index]

vis = stx_calibrate_visibility(vis)

;;******* Compute the Back Projection map
bp_nat_map = stx_bproj(vis,imsize,pixel,aux_data)

Expand Down
47 changes: 47 additions & 0 deletions stix/idl/processing/imaging/stx_uv_points.pro
Original file line number Diff line number Diff line change
@@ -0,0 +1,47 @@
;+
;
; NAME:
; stx_uv_points
;
; PURPOSE:
; Compute the coordinates of the (u,v) points sampled by specific STIX sub-collimators. For more information, we refer to:
;
; - Giordano et al., "The Process of Data Formation for the Spectrometer/Telescope for Imaging X-rays
; (STIX) in Solar Orbiter", SIAM Journal on Imaging Sciences, 2015
;
; - Massa et al., "STIX imaging I - Concept", arxiv, 2023
;
; CALLING SEQUENCE:
; uv_points = stx_uv_points(subc_index)
;
; INPUTS:
; subc_index: array containing the indices of the subcollimators to be considered for computing the (u,v) points
;
; KEYWORDS:
; d_sep: distance between the front and the rear grid [in mm]. Default, 545.30
;
; d_det: distance between the rear grid and the detector [in mm]. Default, 47.78
;
; HISTORY: July 2023, Massa P. (WKU), based on 'stx_uv_points_giordano'
;
; CONTACT:
; paolo.massa@wku.edu
;-

function stx_uv_points, subc_index, d_sep=d_sep, d_det=d_det

default, d_sep, 545.30
default, d_det, 47.78

subc_str = stx_construct_subcollimator()
subc_str = subc_str(subc_index)

factor_f = (d_sep + d_det) / subc_str.front.pitch / 3600.0d / !RADEG ;; in arcsec^-1
factor_r = d_det / subc_str.rear.pitch / 3600.0d / !RADEG ;; in arcsec^-1

u = -subc_str.PHASE * (cos(subc_str.front.angle * !dtor) * factor_f - cos(subc_str.rear.angle * !dtor) * factor_r)
v = -subc_str.PHASE * (sin(subc_str.front.angle * !dtor) * factor_f - sin(subc_str.rear.angle * !dtor) * factor_r)

return, {u:u, v:v}

end
17 changes: 0 additions & 17 deletions stix/idl/processing/imaging/stx_uv_points_giordano.pro

This file was deleted.

21 changes: 6 additions & 15 deletions stix/idl/processing/vis/stx_calibrate_visibility.pro
Original file line number Diff line number Diff line change
Expand Up @@ -19,11 +19,9 @@
;
; KEYWORDS:
;
; phase_calib_factors: 32-element array containing the phase calibration factors for each detectors (degrees).
; phase_calib_factors: 32-element array containing the phase calibration factors for each detector (degrees).
; The default phase calibration factors consist of four terms:
; - a grid correction factor, which keeps into account the phase of the front and the rear grid;
; - a projection correction factor, if the 'xy_flare' estimate of the flare location is provided
; in the visibility structure;
; - an "ad hoc" phase correction factor, which removes systematic residual errors.
; The cause of these systematic errors is still under investigation;
; - a factor which is added so that the reconstructed image is centered in the
Expand All @@ -46,6 +44,8 @@
; Calibrated 'stx_visibility' structure.
;
; HISTORY: August 2022, Massa P., created
; July 2023, Massa P., removed visibility phase 'projection correction' since the new definition of
; (u,v)-points is adopted (see stx_uv_points).
;
; CONTACT:
; paolo.massa@wku.edu
Expand All @@ -63,26 +63,17 @@ modulation_efficiency = !pi^3./(8.*sqrt(2.))

;; Grid phase correction
tmp = read_csv(loc_file( 'GridCorrection.csv', path = getenv('STX_VIS_PHASE') ), header=header, table_header=tableheader, n_table_header=2 )
grid_phase_corr = tmp.field2[vis.ISC - 1]; * (-vis.phase_sense)

;; Projection correction factor
xy_flare = vis[0].XY_FLARE
phase_proj_corr = fltarr(n_vis)
if ~xy_flare[0].isnan() then begin
proj_corr_factor = -xy_flare[0] * 360. * !pi / (180. * 3600. * 8.8) * (r2d_sep + f2r_sep/2.)
phase_proj_corr = phase_proj_corr + proj_corr_factor
;phase_proj_corr = phase_proj_corr * (-vis.phase_sense)
endif
grid_phase_corr = tmp.field2[vis.ISC - 1]

;; "Ad hoc" phase correction (for removing residual errors)
tmp = read_csv(loc_file( 'PhaseCorrFactors.csv', path = getenv('STX_VIS_PHASE')), header=header, table_header=tableheader, n_table_header=3 )
ad_hoc_phase_corr = tmp.field2[vis.ISC - 1]; * (-vis.phase_sense)
ad_hoc_phase_corr = tmp.field2[vis.ISC - 1]

;; Mapcenter correction
phase_mapcenter_corr = -2 * !pi * (vis.XYOFFSET[0] * vis.U + vis.XYOFFSET[1] * vis.V ) * !radeg

default, amp_calib_factors, fltarr(n_vis) + modulation_efficiency
default, phase_calib_factors, grid_phase_corr + ad_hoc_phase_corr + phase_proj_corr + phase_mapcenter_corr
default, phase_calib_factors, grid_phase_corr + ad_hoc_phase_corr + phase_mapcenter_corr
default, syserr_sigamp, 0.05 ;; 5% systematic error in the visibility amplitudes; arbitrary choice

if vis[0].CALIBRATED eq 1 then message, "This visibility structure is already calibrated"
Expand Down
14 changes: 4 additions & 10 deletions stix/idl/processing/vis/stx_pixel_data_summed2visibility.pro
Original file line number Diff line number Diff line change
Expand Up @@ -52,6 +52,7 @@
; - CALIBRATED: 0 if the values of the visibility amplitudes and phases are not calibrated, 1 otherwise
;
; HISTORY: August 2022, Massa P., created
; July 2023, Massa P., made it compatible with the new definition of (u,v)-points (see stx_uv_points)
;
; CONTACT:
; paolo.massa@wku.edu
Expand All @@ -77,16 +78,9 @@ subc_str = subc_str[subc_index]

;;************** Define (u,v) points

; take average of front and rear grid pitches (mm)
pitch = (subc_str.front.pitch + subc_str.rear.pitch) / 2.0d
; take average of front and rear grid orientation
orientation = (subc_str.front.angle + subc_str.rear.angle) / 2.0d
; convert pitch from mm to arcsec
pitch = pitch / f2r_sep * 3600.0d * !RADEG
; calculate u and v
uv = 1.0 / pitch
u = uv * cos(orientation * !DTOR) * (-subc_str.PHASE) ; TO BE REMOVED!
v = uv * sin(orientation * !DTOR) * (-subc_str.PHASE) ; TO BE REMOVED!
uv_data = stx_uv_points(subc_index)
u = uv_data.u
v = uv_data.v

;;************** Define visibility values

Expand Down