-
Notifications
You must be signed in to change notification settings - Fork 3
/
nh3_testcubes.py
218 lines (201 loc) · 9.26 KB
/
nh3_testcubes.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
import pyspeckit.spectrum.models.ammonia as ammonia
import pyspeckit.spectrum.models.ammonia_constants as nh3con
from pyspeckit.spectrum.units import SpectroscopicAxis as spaxis
import os
import sys
import numpy as np
import astropy.units as u
from astropy.io import fits
from spectral_cube import SpectralCube
from astropy.utils.console import ProgressBar
from astropy import log
log.setLevel('ERROR')
def generate_cubes(nCubes=100, nBorder=1, noise_rms=0.1,
output_dir='random_cubes', fix_vlsr=True,
random_seed=None):
"""
This places nCubes random cubes into the specified output directory
"""
if not os.path.isdir(output_dir):
os.mkdir(output_dir)
xarr11 = spaxis((np.linspace(-500, 499, 1000) * 5.72e-6
+ nh3con.freq_dict['oneone'] / 1e9),
unit='GHz',
refX=nh3con.freq_dict['oneone'] / 1e9,
velocity_convention='radio', refX_unit='GHz')
xarr22 = spaxis((np.linspace(-500, 499, 1000) * 5.72e-6
+ nh3con.freq_dict['twotwo'] / 1e9), unit='GHz',
refX=nh3con.freq_dict['twotwo'] / 1e9,
velocity_convention='radio', refX_unit='GHz')
nDigits = int(np.ceil(np.log10(nCubes)))
if random_seed:
np.random.seed(random_seed)
nComps = np.random.choice([1, 2], nCubes)
Temp1 = 8 + np.random.rand(nCubes) * 17
Temp2 = 8 + np.random.rand(nCubes) * 17
if fix_vlsr:
Voff1 = np.zeros(nCubes)
else:
Voff1 = np.random.rand(nCubes) * 5 - 2.5
Voff2 = Voff1 + np.random.rand(nCubes) * 5 - 2.5
logN1 = 13 + 1.5 * np.random.rand(nCubes)
logN2 = 13 + 1.5 * np.random.rand(nCubes)
Width1NT = 0.1 * np.exp(1.5 * np.random.randn(nCubes))
Width2NT = 0.1 * np.exp(1.5 * np.random.randn(nCubes))
# Width1 = 0.08 + 1.0 * np.random.rand(nCubes)
# Width2 = 0.08 + 1.0 * np.random.rand(nCubes)
Width1 = np.sqrt(Width1NT + 0.08**2)
Width2 = np.sqrt(Width2NT + 0.08**2)
# Find where centroids are too close
too_close = np.where(np.abs(Voff1-Voff2)<np.max(np.column_stack((Width1, Width2)), axis=1))
# Move the centroids farther apart by the length of largest line width
min_Voff = np.min(np.column_stack((Voff2[too_close],Voff1[too_close])), axis=1)
max_Voff = np.max(np.column_stack((Voff2[too_close],Voff1[too_close])), axis=1)
Voff1[too_close]=min_Voff-np.max(np.column_stack((Width1[too_close], Width2[too_close])), axis=1)/2.
Voff2[too_close]=max_Voff+np.max(np.column_stack((Width1[too_close], Width2[too_close])), axis=1)/2.
scale = np.array([[0.2, 0.1, 0.5, 0.01]])
gradX1 = np.random.randn(nCubes, 4) * scale
gradY1 = np.random.randn(nCubes, 4) * scale
gradX2 = np.random.randn(nCubes, 4) * scale
gradY2 = np.random.randn(nCubes, 4) * scale
params1 = [{'ntot':14,
'width':1,
'xoff_v':0.0}] * nCubes
params2 = [{'ntot':14,
'width':1,
'xoff_v':0.0}] * nCubes
hdrkwds = {'BUNIT': 'K',
'INSTRUME': 'KFPA ',
'BMAJ': 0.008554169991270138,
'BMIN': 0.008554169991270138,
'TELESCOP': 'GBT',
'WCSAXES': 3,
'CRPIX1': 2,
'CRPIX2': 2,
'CRPIX3': 501,
'CDELT1': -0.008554169991270138,
'CDELT2': 0.008554169991270138,
'CDELT3': 5720.0,
'CUNIT1': 'deg',
'CUNIT2': 'deg',
'CUNIT3': 'Hz',
'CTYPE1': 'RA---TAN',
'CTYPE2': 'DEC--TAN',
'CTYPE3': 'FREQ',
'CRVAL1': 0.0,
'CRVAL2': 0.0,
'LONPOLE': 180.0,
'LATPOLE': 0.0,
'EQUINOX': 2000.0,
'SPECSYS': 'LSRK',
'RADESYS': 'FK5',
'SSYSOBS': 'TOPOCENT'}
truekwds = ['NCOMP', 'LOGN1', 'LOGN2', 'VLSR1', 'VLSR2',
'SIG1', 'SIG2', 'TKIN1', 'TKIN2']
for i in ProgressBar(range(nCubes)):
xmat, ymat = np.indices((2 * nBorder + 1, 2 * nBorder + 1))
cube11 = np.zeros((xarr11.shape[0], 2 * nBorder + 1, 2 * nBorder + 1))
cube22 = np.zeros((xarr22.shape[0], 2 * nBorder + 1, 2 * nBorder + 1))
Tmax11a, Tmax11b, Tmax22a, Tmax22b = (0,) * 4
for xx, yy in zip(xmat.flatten(), ymat.flatten()):
T1 = Temp1[i] * (1 + gradX1[i, 0] * (xx - 1)
+ gradY1[i, 0] * (yy - 1)) + 5
T2 = Temp2[i] * (1 + gradX2[i, 0] * (xx - 1)
+ gradY2[i, 0] * (yy - 1)) + 5
if T1 < 2.74:
T1 = 2.74
if T2 < 2.74:
T2 = 2.74
W1 = np.abs(Width1[i] * (1 + gradX1[i, 1] * (xx - 1)
+ gradY1[i, 1] * (yy - 1)))
W2 = np.abs(Width2[i] * (1 + gradX2[i, 1] * (xx - 1)
+ gradY2[i, 1] * (yy - 1)))
V1 = Voff1[i] + (gradX1[i, 2] * (xx - 1) + gradY1[i, 2] * (yy - 1))
V2 = Voff2[i] + (gradX2[i, 2] * (xx - 1) + gradY2[i, 2] * (yy - 1))
N1 = logN1[i] * (1 + gradX1[i, 3] * (xx - 1)
+ gradY1[i, 3] * (yy - 1))
N2 = logN2[i] * (1 + gradX2[i, 3] * (xx - 1)
+ gradY2[i, 3] * (yy - 1))
if nComps[i] == 1:
spec11 = ammonia.cold_ammonia(xarr11, T1,
ntot=N1,
width=W1,
xoff_v=V1)
spec22 = ammonia.cold_ammonia(xarr22, T1,
ntot=N1,
width=W1,
xoff_v=V1)
if (xx == nBorder) and (yy == nBorder):
Tmax11a = np.max(spec11)
Tmax22a = np.max(spec22)
Tmax11 = np.max(spec11)
Tmax22 = np.max(spec22)
if nComps[i] == 2:
spec11a = ammonia.cold_ammonia(xarr11, T1,
ntot=N1,
width=W1,
xoff_v=V1)
spec11b = ammonia.cold_ammonia(xarr11, T2,
ntot=N2,
width=W2,
xoff_v=V2)
spec11 = spec11a + spec11b
spec22a = ammonia.cold_ammonia(xarr22, T1,
ntot=N1,
width=W1,
xoff_v=V1)
spec22b = ammonia.cold_ammonia(xarr22, T2,
ntot=N2,
width=W2,
xoff_v=V2)
spec22 = spec22a + spec22b
if (xx == nBorder) and (yy == nBorder):
Tmax11a = np.max(spec11a)
Tmax11b = np.max(spec11b)
Tmax22a = np.max(spec22a)
Tmax22b = np.max(spec22b)
Tmax11 = np.max(spec11)
Tmax22 = np.max(spec22)
cube11[:, yy, xx] = spec11
cube22[:, yy, xx] = spec22
cube11 += np.random.randn(*cube11.shape) * noise_rms
cube22 += np.random.randn(*cube22.shape) * noise_rms
hdu11 = fits.PrimaryHDU(cube11)
for kk in hdrkwds:
hdu11.header[kk] = hdrkwds[kk]
for kk, vv in zip(truekwds, [nComps[i], logN1[i], logN2[i],
Voff1[i], Voff2[i], Width1[i], Width2[i],
Temp1[i], Temp2[i]]):
hdu11.header[kk] = vv
hdu11.header['TMAX'] = Tmax11
hdu11.header['TMAX-1'] = Tmax11a
hdu11.header['TMAX-2'] = Tmax11b
hdu11.header['RMS'] = noise_rms
hdu11.header['CRVAL3'] = 23694495500.0
hdu11.header['RESTFRQ'] = 23694495500.0
hdu11.writeto(output_dir + '/random_cube_NH3_11_'
+ '{0}'.format(i).zfill(nDigits)
+ '.fits',
overwrite=True)
hdu22 = fits.PrimaryHDU(cube22)
for kk in hdrkwds:
hdu22.header[kk] = hdrkwds[kk]
for kk, vv in zip(truekwds, [nComps[i], logN1[i], logN2[i],
Voff1[i], Voff2[i], Width1[i], Width2[i],
Temp1[i], Temp2[i]]):
hdu22.header[kk] = vv
hdu22.header['TMAX'] = Tmax22
hdu22.header['TMAX-1'] = Tmax22a
hdu22.header['TMAX-2'] = Tmax22b
hdu22.header['RMS'] = noise_rms
hdu22.header['CRVAL3'] = 23722633600.0
hdu22.header['RESTFRQ'] = 23722633600.0
hdu22.writeto(output_dir + '/random_cube_NH3_22_'
+ '{0}'.format(i).zfill(nDigits) + '.fits',
overwrite=True)
if __name__ == '__main__':
print(sys.argv)
if len(sys.argv) > 1:
generate_cubes(nCubes=int(sys.argv[1]))
else:
generate_cubes()