pro rd_specline_hitran, fp, ln ;+ ; NAME: ; rd_specline_hitran ; PURPOSE: (one line) ; Read the next line from an open HITRAN file, and fill a spectral line structure ; DESCRIPTION: ; CATEGORY: ; Spectra ; CALLING SEQUENCE: ; rd_specline_hitran, file_pointer, spectral_line ; INPUTS: ; file_pointer - a file pointer to an open HITRAN file. ; spectral_line - a previously defined spectral line structure ; OPTIONAL INPUT PARAMETERS: ; none ; KEYWORD INPUT PARAMETERS: ; none ; KEYWORD OUTPUT PARAMETERS: ; none ; OUTPUTS: ; Fills the spectral line, spectral_line ; COMMON BLOCKS: ; None ; SIDE EFFECTS: ; RESTRICTIONS: ; spectral_line must be defined before calling. ; PROCEDURE: ; MODIFICATION HISTORY: ; Written 2000 October, by Leslie Young, SwRI ; 2008 Dec 16 LAY Modified for HITRAN2004 ; 2013 Aug 08 LAY Modified for HITRAN2012 ;- ON_ERROR,2 ; Return to caller if an error occurs ln.mol=-1 ln.nu0 = 0. fs=fstat(fp) if fs.open eq 0 or fs.read eq 0 then return if eof(fp) then return case !spectralib.version of 'HITRAN96': format='(I2,I1,F12.6,2E10.3,2F5.4,F10.4,F4.2,F8.6,2A3, 2A9, 3I1,3I2)' 'HITRAN2004': format='(I2,I1,F12.6,2E10.3,2F5.4,F10.4,F4.2,F8.6,2A15,2A15,6I1,6I2'$ + ',A1, 2F7.1)' 'HITRAN2008': format='(I2,I1,F12.6,2E10.3,2F5.4,F10.4,F4.2,F8.6,2A15,2A15,6I1,6I2'$ + ',A1, 2F7.1)' 'HITRAN2012': format='(I2,I1,F12.6,2E10.3,2F5.4,F10.4,F4.2,F8.6,2A15,2A15,6I1,6I2'$ + ',A1, 2F7.1)' else: format='(I2,I1,F12.6,2E10.3,2F5.4,F10.4,F4.2,F8.6,2A3,2A9,3I1,3I2)' endcase junk='' mol=0 iso=0 nu0=0.0d s=0.0d a=0.0d r2=0.0d gamma_air=0.0d gamma_self=0.0d E_lower=0.0d n=0.0d delta=0.0d case !spectralib.version of 'HITRAN96': v_def = 0 'HITRAN2004': v_def = ' ' 'HITRAN2008': v_def = ' ' 'HITRAN2012': v_def = ' ' else: v_def = ' ' endcase v_upper=v_def v_lower=v_def Q_upper=' ' Q_lower=' ' ier0=0 ier1=0 ier2=0 ier3=0 ier4=0 ier5=0 iref0=0 iref1=0 iref2=0 iref3=0 iref4=0 iref5=0 linemix_flag = ' ' g_upper = 0.0d g_lower = 0.0d line='' readf,fp,junk if strlen(junk) le 99 then readf,fp,junk case !spectralib.version of 'HITRAN96': $ reads, junk, $ mol, iso, nu0, s, r2, gamma_air, gamma_self, E_lower, n, delta, v_upper, v_lower, $ Q_upper, Q_lower, ier0,ier1,ier2, iref0,iref1,iref2, format=format 'HITRAN2004': $ reads, junk, $ mol, iso, nu0, s, a, gamma_air, gamma_self, E_lower, n, delta, v_upper, v_lower, $ Q_upper, Q_lower, ier0,ier1,ier2,ier3,ier4,ier5, $ iref0,iref1,iref2,iref3,iref4,iref5,linemix_flag, g_upper, g_lower, $ format=format 'HITRAN2008': $ reads, junk, $ mol, iso, nu0, s, a, gamma_air, gamma_self, E_lower, n, delta, v_upper, v_lower, $ Q_upper, Q_lower, ier0,ier1,ier2,ier3,ier4,ier5, $ iref0,iref1,iref2,iref3,iref4,iref5,linemix_flag, g_upper, g_lower, $ format=format 'HITRAN2012': $ reads, junk, $ mol, iso, nu0, s, a, gamma_air, gamma_self, E_lower, n, delta, v_upper, v_lower, $ Q_upper, Q_lower, ier0,ier1,ier2,ier3,ier4,ier5, $ iref0,iref1,iref2,iref3,iref4,iref5,linemix_flag, g_upper, g_lower, $ format=format else: endcase ln.mol=mol ln.iso=iso ln.nu0=nu0 ln.s=s ln.a = a ln.r2 = r2 ln.gamma_air=gamma_air ln.gamma_self=gamma_self ln.E_lower=E_lower ln.n=n ln.delta=delta ln.v_upper=v_upper ln.v_lower=v_lower ln.Q_upper=Q_upper ln.Q_lower=Q_lower ln.ier=[ier0,ier1,ier2,ier3,ier4,ier5] ln.iref=[iref0,iref1,iref2,iref3,iref4,iref5] ln.linemix_flag = linemix_flag ln.g_upper = g_upper ln.g_lower = g_lower ln.t=296. ln.t_orig=ln.t ln.p = 1.01325e6 ln.dopw = dopwidth(nu0, ln.t, molec_weight(mol,iso)) ln.p_orig = ln.p ln.s_orig=s ln.gamma_air_orig=gamma_air ln.gamma_self_orig=gamma_self end