Документ взят из кэша поисковой машины. Адрес оригинального документа : http://star.arm.ac.uk/~csj/software_store/idlines/idlines.pro
Дата изменения: Mon Jun 2 19:44:55 2003
Дата индексирования: Tue Oct 2 10:27:11 2012
Кодировка:

Поисковые слова: п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п п
;---------------------------------------------------------------------------
;
; IDL procedures for plotting spectra and line identifiers
;
;
; idlines ; example of how to use tools
;
; readspec : read spectrum from 2 column x-y file with 3 line header
; readspecfit : read spectrum and fit from 4 column file with " " "
;
; plotspec : plot a spectrum onto a new frame
; [plotfit : plot a spectrum onto the current frame - redundant]
;
; readlist : read linelist in lte_lines format
;
; marklist : marks all lines in linelist
; markelem : marks positions, identities and wavelengths of element
; markspec : marks positions of spectral lines from a given element
;
;
;---------------------------------------------------------------------------
;
; For additional routines, see also
;
; see /home/csj/idl/...
; /home/csj/idl/dipsotools.pro
;
;---------------------------------------------------------------------------
;
; Author: C S Jeffery, Armagh Observatory, 2000.
;
;---------------------------------------------------------------------------




PRO idlines, specfile, linelist, xlo, xhi, atoms, fit

;--------------------------------------------------------------------------
;
; idlines
;
; IDL proecudure to diplsay a region of spectrum, together with
; a fit, at full scale and 10 x full scale. Positions of lines in
; common atomic spectra up to argon are marked in the lower
; (full-scale) panel. Lines of sleceted atoms are identified in the
; middle panel, the enlarged spectrum is shown in the uppermost panel.
;
; Command arguments:
;
; specfile - character string - name of file containing spectrum (and fit)
; assumed to be x,y file with 3-line header
; linelist - character string - name of file containing linelist
; assumed to be in 'lte_lines' format
; xlo - real variable - lower limit of plot
; xhi - real variable - upper limit of plot
; atoms - integer vector - atomic numbers of spectra to be identified
; fit - integer variable - switch = 0 don't read or plot fit
; = 1 read and plot fit from specfile
;
;---------------------------------------------------------------------------


;; get the spectrum to plot

IF fit EQ 1 THEN BEGIN
readspecfit, specfile, specnw, specwav, specflx, specfit
ENDIF ELSE BEGIN
readspec, specfile, specnw, specwav, specflx
ENDELSE

;; get the linelist with identifications

readlist, linelist, listn, listwav, listlab, listion, listmlt, $
listgf, listep



;; plot the spectrum and fit

plotspec, specwav, specflx, xlo, xhi, 0.0, 2.1

IF fit EQ 1 THEN $
oplot, specwav, specfit, color=127 ; overplot the fit

;; mark all lines
;; marklist, xlo, xhi, listn, listwav, listion, listlab, listmlt


;; mark individual lines

FOR zel = 0,n_elements(atoms)-1 DO BEGIN
markelem, atoms(zel), xlo, xhi, listn, listwav, listion, listlab, listmlt
ENDFOR


;; mark separate spectra

FOR zel = 2,26 DO BEGIN
markspec, zel, xlo, xhi, listn, listwav, listion, listlab, listmlt, $
listgf
ENDFOR



;; expand the plot

oplot, [xlo, xhi], [1.5,1.5]
oplot, [xlo, xhi], [2.0,2.0], linestyle=1
oplot, specwav, 10 * specflx - 8, clip = [xlo,1.5,xhi,2.1]
IF fit EQ 1 THEN $
oplot, specwav, 10 * specfit - 8, clip = [xlo,1.5,xhi,2.1], color=127


END






PRO plotspec, specwav, specflx, xlo, xhi, ylo, yhi

; establish the box

plot, [xlo,xhi],[ylo, yhi], $
xrange=[xlo,xhi], yrange=[ylo,yhi], xstyle = 1, ystyle=1, /nodata

oplot, specwav, specflx ; plot the spectrum
oplot, [xlo,xhi], [1,1], linestyle = 1 ; mark the continuum

END




PRO plotfit, specwav, specfit, color

!p.color = color
oplot, specwav, specfit ; overplot the fit

END






PRO markelem, zel, xlo, xhi, listn, listwav, listion, listlab, listmlt

FOR i = 0,listn-1 DO BEGIN
IF zel EQ FIX (listion(i)) THEN BEGIN
wl = listwav(i)
IF xlo LT wl AND wl LT xhi THEN BEGIN


oplot, [wl,wl], [1.05,1.15],linestyle = 0,color=159 ; mark line position
oplot, [wl,wl], [1.45,1.55],linestyle = 0,color=159 ; mark line position

lab = listlab(i)
xyouts, wl,1.18, lab, orientation=90, noclip=0, color=223, charsize=0.7

zmu = listmlt(i)
xmu = fix(zmu)
IF zmu-xmu LT 0.001 THEN wmu = string(xmu,format='(i3)') $
ELSE wmu = string(zmu,format='(f6.2)')
xyouts, wl,1.25, wmu , orientation=90, noclip=0, color=223, charsize=0.5

wls = string(listwav(i),format='(f7.2)')
xyouts, wl,1.32, wls , orientation=90, noclip=0, color=223, charsize=0.5

ENDIF
ENDIF
ENDFOR

END




PRO markspec, zel, xlo, xhi, listn, listwav, listion, listlab, listmlt, listgf

; Data for constructing labels

elem = StrArr(50)
elem( 1: 2) = [ 'H ', 'He' ]
elem( 3:10) = [ 'Be', 'Li', 'B ', 'C ', 'N ', 'O ', 'F ', 'Ne' ]
elem(11:18) = [ 'Na', 'Mg', 'Al', 'Si', 'P ', 'S ', 'Cl', 'A ' ]
elem(19:20) = [ 'K ', 'Ca' ]
elem(31:36) = [ 'Ga', 'Ge', 'As', 'Se', 'Br', 'Kr' ]
elem(21:30) = [ 'Sc', 'Ti', 'V ', 'Cr', 'Mn', 'Fe', 'Co', 'Ni', 'Cu', 'Zn' ]

yzel = FltArr(50) ; position of each atom
yzel( 1: 2) = [ 0.00, 0.02 ]
yzel( 3:10) = [ 0.04, 0.05, 0.06, 0.08, 0.12, 0.16, 0.18, 0.20 ]
yzel(11:18) = [ 0.22, 0.24, 0.26, 0.28, 0.30, 0.32, 0.34, 0.36 ]
yzel(19:20) = [ 0.38, 0.40 ]
yzel(31:36) = [ 1.04, 1.06, 1.08, 1.10, 1.12, 1.14 ]
yzel(21:30) = [ 0.42, 0.44, 0.48, 0.50, 0.52, 0.56, 0.58, 0.60, 0.62, 0.64 ]

lzel = IntArr(50) ; linestyles for each atom
lzel( 1: 2) = [ 0, 0 ]
lzel( 3:10) = [ 1, 1, 1, 0, 0, 0, 1, 0 ]
lzel(11:18) = [ 1, 0, 1, 0, 1, 0, 1, 0 ]
lzel(19:20) = [ 1, 0 ]
lzel(31:36) = [ 1, 1, 1, 1, 1, 1 ]
lzel(21:30) = [ 1, 0, 1, 0, 1, 0, 1, 0, 1, 1 ]

; level at which this ion is to be marked
yb = yzel(zel)
nid = 0

gfmax = -99.
FOR i = 0,listn-1 DO BEGIN
IF zel EQ FIX (listion(i)) THEN BEGIN
wl = listwav(i)
IF xlo LT wl AND wl LT xhi THEN BEGIN
gfmax = max([gfmax,listgf(i)])
ENDIF
ENDIF
ENDFOR
gfmax = alog10(gfmax)

FOR i = 0,listn-1 DO BEGIN
IF zel EQ FIX (listion(i)) THEN BEGIN
wl = listwav(i)
IF xlo LT wl AND wl LT xhi THEN BEGIN
dy = 0.03
dy = 0.025 / (gfmax - alog10(listgf(i)) + 0.5) + 0.001
oplot, [wl,wl], [yb,yb+dy], linestyle = lzel(zel), color=159
nid = nid + 1
ENDIF
ENDIF
ENDFOR

IF nid GT 0 THEN BEGIN
oplot, [xlo,xhi], [yb,yb], linestyle = 0,color=159
dx = (xhi-xlo)/20.
xyouts, xhi-(1+lzel(zel))*dx , yb+0.005, elem(zel), noclip=0, color=191
ENDIF

END





;----------------------------------------------------------------------------
;
; read lines from a linelist and returns them as an array of wavelengths,
; labels, and ion identifications
;
;----------------------------------------------------------------------------

PRO readlist, file, listn, listwav, listlab, listion, listmlt, listgf, listep


; Data for constructing labels

elem = StrArr(50)
elem( 1: 2) = [ 'H ', 'He' ]
elem( 3:10) = [ 'Be', 'Li', 'B ', 'C ', 'N ', 'O ', 'F ', 'Ne' ]
elem(11:18) = [ 'Na', 'Mg', 'Al', 'Si', 'P ', 'S ', 'Cl', 'A ' ]
elem(19:20) = [ 'K ', 'Ca' ]
elem(31:36) = [ 'Ga', 'Ge', 'As', 'Se', 'Br', 'Kr' ]
elem(21:30) = [ 'Sc', 'Ti', 'V ', 'Cr', 'Mn', 'Fe', 'Co', 'Ni', 'Cu', 'Zn' ]

ion = StrArr(6)
ion = [ ' I', ' II', 'III', ' IV', ' V', ' VI' ]

; Data array definitions

nmax = 32000
listwav = FltArr(nmax)
listlab = StrArr(nmax)
listion = FltArr(nmax)
listmlt = FltArr(nmax)
listgf = FltArr(nmax)
listep = FltArr(nmax)

; ---


; open the data file

OPENR, lulist, file, /get_lun


; skip over atomic data

element = ''
WHILE StrPos(element,'ENDATOMS') LT 0 DO BEGIN
READF, lulist, element
ENDWHILE

; read individual lines in list

i = 0 & valid = 0 & nz = 999
WHILE nz NE -1 DO BEGIN
ON_IOERROR, end_of_data
READF, lulist, nz, ni, wl, gf, ed, rd, ep, mult
IF nz GT 0 THEN BEGIN
listwav(i) = wl
listion(i) = nz + (ni-1)/100.
listlab(i) = elem(nz) + ion(ni-1)
listmlt(i) = mult
listgf(i) = gf
listep(i) = ep
; print, i, ' ', listlab(i), listmlt(i), listwav(i), listmlt(i)
i = i + 1
ENDIF
ENDWHILE

; end of data, close file

end_of_data:
CLOSE, lulist
FREE_LUN, lulist


; tidy up output arrays

listn = i
listwav = listwav(0:i-1)
listion = listion(0:i-1)
listlab = listlab(0:i-1)
listmlt = listmlt(0:i-1)
listgf = listgf(0:i-1)
listep = listep(0:i-1)

END





;----------------------------------------------------------------------------
;
; read spectrum and fit (sfit_synth output)
;
;----------------------------------------------------------------------------

PRO readspecfit, file, specnw, specwav, specflx, specfit


; Data array definitions

smax = 32000
specwav = FltArr(smax)
specflx = FltArr(smax)
specfit = FltArr(smax)

; ---


; open the data file

OPENR, luspec, file, /get_lun


; skip over header

header = ''
FOR i = 1,3 DO BEGIN
READF, luspec, header
print, header
ENDFOR

; read individual lines in list

i = 0 & valid = 0
WHILE valid EQ 0 DO BEGIN
ON_IOERROR, end_of_data
READF, luspec, wl, fl, fs, ff
specwav(i) = wl
specflx(i) = fl
specfit(i) = ff
i = i + 1
ENDWHILE

; end of data, close file

end_of_data:
CLOSE, luspec
FREE_LUN, luspec


; tidy up output arrays

specnw = i
specwav = specwav(0:i-1)
specflx = specflx(0:i-1)
specfit = specfit(0:i-1)

END



;----------------------------------------------------------------------------
;
; read spectrum only (free format)
;
;----------------------------------------------------------------------------

PRO readspec, file, specnw, specwav, specflx


; Data array definitions

specwav = FltArr(32000)
specflx = FltArr(32000)

; ---


; open the data file

OPENR, luspec, file, /get_lun


; skip over header

header = ''
FOR i = 1,3 DO BEGIN
READF, luspec, header
ENDFOR

; read individual lines in list

i = 0 & valid = 0
WHILE valid EQ 0 DO BEGIN
ON_IOERROR, end_of_data
READF, luspec, wl, fl
specwav(i) = wl
specflx(i) = fl
i = i + 1
ENDWHILE

; end of data, close file

end_of_data:
CLOSE, luspec
FREE_LUN, luspec


; tidy up output arrays

specnw = i
specwav = specwav(0:i-1)
specflx = specflx(0:i-1)

END