Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
ift
nifty_gridder
Commits
ae8d400c
Commit
ae8d400c
authored
Jun 05, 2019
by
Martin Reinecke
Browse files
switch to pixel size
parent
6e72becc
Changes
2
Hide whitespace changes
Inline
Side-by-side
demo.py
View file @
ae8d400c
...
...
@@ -27,8 +27,8 @@ import nifty_gridder as ng
# dirty: float(nxdirty, nydirty): the dirty image
f0
=
1e9
# rough observation frequency
fov
=
np
.
pi
/
180
/
60
# assume 1 arcmin FOV
npixdirty
=
1024
pixsize
=
np
.
pi
/
180
/
60
/
npixdirty
# assume 1 arcmin FOV
speedoflight
=
3e8
# number of rows in the measurement set
...
...
@@ -38,23 +38,24 @@ nrow = 10000
nchan
=
100
# Frequency for all channels in the measurement set [Hz].
freq
=
f0
+
np
.
arange
(
nchan
)
*
(
f0
/
nchan
)
# just use silly values for this example
# Invent mock UVW data. For every row, there is one UVW triple
# NOTE: I don't know how to set the w values properly, so I use zeros.
uvw
=
(
np
.
random
.
rand
(
nrow
,
3
)
-
0.5
)
/
(
fov
/
npixdirty
*
f0
/
speedoflight
)
uvmax
=
pixsize
*
np
.
max
(
freq
)
/
sp
uvw
=
(
np
.
random
.
rand
(
nrow
,
3
)
-
0.5
)
/
(
pixsize
*
f0
/
speedoflight
)
uvw
[:,
2
]
=
0.
# Frequency for all channels in the measurement set [Hz].
freq
=
f0
+
np
.
arange
(
nchan
)
*
(
f0
/
nchan
)
# just use silly values for this example
# Build Baselines object from the geometrical information
baselines
=
ng
.
Baselines
(
coord
=
uvw
,
freq
=
freq
)
# Build GridderConfig object describing how gridding should be done
#
fov_x and fov
_y are given in radians
gconf
=
ng
.
GridderConfig
(
nxdirty
=
npixdirty
,
nydirty
=
npixdirty
,
epsilon
=
1e-7
,
fov_x
=
fov
,
fov_y
=
fov
)
#
pixsize_x and pixsize
_y are given in radians
gconf
=
ng
.
GridderConfig
(
nxdirty
=
npixdirty
,
nydirty
=
npixdirty
,
epsilon
=
1e-7
,
pixsize_x
=
pixsize
,
pixsize_y
=
pixsize
)
# At this point everything about the experimental setup is known and explained
...
...
nifty_gridder.cc
View file @
ae8d400c
...
...
@@ -393,9 +393,9 @@ template<typename T> class GridderConfig
public:
GridderConfig
(
size_t
nxdirty
,
size_t
nydirty
,
double
epsilon
,
double
fov
_x
,
double
fov
_y
)
double
pixsize
_x
,
double
pixsize
_y
)
:
nx_dirty
(
nxdirty
),
ny_dirty
(
nydirty
),
ucorr
(
fov_x
/
nx_dirty
),
vcorr
(
fov_y
/
ny_dirt
y
),
ucorr
(
pixsize_x
),
vcorr
(
pixsize_
y
),
w
(
get_w
(
epsilon
)),
nsafe
((
w
+
1
)
/
2
),
nu
(
max
(
2
*
nsafe
,
2
*
nx_dirty
)),
nv
(
max
(
2
*
nsafe
,
2
*
ny_dirty
)),
beta
(
2.3
*
w
),
...
...
@@ -406,8 +406,8 @@ template<typename T> class GridderConfig
myassert
((
nx_dirty
&
1
)
==
0
,
"nx_dirty must be even"
);
myassert
((
ny_dirty
&
1
)
==
0
,
"ny_dirty must be even"
);
myassert
(
epsilon
>
0
,
"epsilon must be positive"
);
myassert
(
fov_x
>
0
,
"fov
_x must be positive"
);
myassert
(
fov_y
>
0
,
"fov
_y must be positive"
);
myassert
(
pixsize_x
>
0
,
"pixsize
_x must be positive"
);
myassert
(
pixsize_y
>
0
,
"pixsize
_y must be positive"
);
auto
tmp
=
correction_factors
(
nu
,
nx_dirty
/
2
+
1
,
w
);
cfu
[
nx_dirty
/
2
]
=
tmp
[
0
];
...
...
@@ -948,10 +948,10 @@ nydirty: int
epsilon: float
required accuracy for the gridding/degridding step
Must be >= 2e-13.
fov
_x: float
Field of view
in x direction (radians)
fov
_y: float
Field of view
in y direction (radians)
pixsize
_x: float
Pixel size
in x direction (radians)
pixsize
_y: float
Pixel size
in y direction (radians)
)"""
;
const
char
*
grid2dirty_DS
=
R"""(
...
...
@@ -1068,7 +1068,7 @@ PYBIND11_MODULE(nifty_gridder, m)
.
def
(
"vis2ms"
,
&
Baselines
<
double
>::
vis2ms
<
complex
<
double
>>
,
BL_vis2ms_DS
,
"vis"
_a
,
"idx"
_a
,
"ms_in"
_a
=
py
::
none
());
py
::
class_
<
GridderConfig
<
double
>>
(
m
,
"GridderConfig"
,
GridderConfig_DS
)
.
def
(
py
::
init
<
size_t
,
size_t
,
double
,
double
,
double
>
(),
"nxdirty"
_a
,
"nydirty"
_a
,
"epsilon"
_a
,
"
fov
_x"
_a
,
"
fov
_y"
_a
)
"nydirty"
_a
,
"epsilon"
_a
,
"
pixsize
_x"
_a
,
"
pixsize
_y"
_a
)
.
def
(
"Nu"
,
&
GridderConfig
<
double
>::
Nu
)
.
def
(
"Nv"
,
&
GridderConfig
<
double
>::
Nv
)
.
def
(
"grid2dirty"
,
&
GridderConfig
<
double
>::
grid2dirty
,
grid2dirty_DS
,
"grid"
_a
)
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment