Scrambler  1
emissions.f90
Go to the documentation of this file.
00001 ! Module for calculating emission of a given transition for a given species
00002 
00003 MODULE Emissions
00004    USE GlobalDeclarations   
00005    IMPLICIT NONE
00006    PUBLIC CalcEmiss
00007 
00008    INTEGER, PUBLIC, PARAMETER :: iOI = 1, iNII = 2, iSII_6716 = 3, iSII_6731 = 4, iHalpha = 5
00009    INTEGER, PUBLIC, PARAMETER :: iSII = 99
00010 
00011 CONTAINS
00012 
00013    ! Main control routine which calculates the emission for given species and transition
00014    SUBROUTINE CalcEmiss(ne,T0,x,nH,iSpecies,emiss)
00015       REAL(KIND=qPREC), INTENT(OUT) :: emiss
00016       REAL(KIND=qPREC), INTENT(IN) :: ne, x, T0, nH   ! electron dens, ion frac, temp, total H dens
00017       INTEGER, INTENT(IN) :: iSpecies                 ! species index
00018 
00019       ! Dimension of b matrix is 5x5 since there are 5 energy levels
00020       REAL(KIND=qPREC) :: b(5,5), g(5), e(5), a(5,5), o(5,5), den(5), indx(5)
00021       REAL(KIND=qPREC) :: levpops, T, emiss1
00022       INTEGER :: lev1, lev2               ! the initial, final transition levels
00023 
00024       ! Use separate routine for H-alpha emission
00025       IF(iSpecies == iHalpha) THEN
00026          CALL alpha_emiss(ne,T0,x,nH,emiss1)
00027          emiss = emiss1
00028          RETURN
00029       END IF 
00030 
00031       ! Put T in 10^4 K
00032       T = T0/(10d0**4d0)
00033 
00034       IF(ne /= 0.) THEN
00035 
00036          SELECT CASE(iSpecies)
00037             CASE(iOI)                    ! [O I]  6300
00038                CALL OI(ne,T,g,e,a,o)
00039                lev1 = 4 ; lev2 = 1
00040             CASE(iNII)                   ! [N II] 6583
00041                CALL NII(ne,T,g,e,a,o)
00042                lev1 = 4 ; lev2 = 3
00043             CASE(iSII_6716)              ! [S II] 6716
00044                CALL SII(ne,T,g,e,a,o)
00045                lev1 = 3 ; lev2 = 1
00046             CASE(iSII_6731)              ! [S II] 6731
00047                CALL SII(ne,T,g,e,a,o)
00048                lev1 = 2 ; lev2 = 1
00049          END SELECT
00050 
00051          CALL Loadb(b,T,ne,g,e,a,o)
00052 
00053          ! Do matrix inversion (values in b and den will be overwritten)
00054          ! Input: matrix b and r.h.s. is (0,0,0,0,1) stored in den
00055          ! Output: LU decomposed matrix in b and density of levels in den
00056          ! From numerical recipes in FORTRAN
00057          den = (/ 0, 0, 0, 0, 1 /)
00058          CALL ludcmp(b,indx)
00059          CALL lubksb(b,indx,den)
00060 
00061       ELSE
00062          den = (/ 1, 0, 0, 0, 0 /)
00063       END IF
00064 
00065       levpops = den(lev1)
00066       IF(levpops == -0.) levpops = 0d0   
00067 
00068       CALL emissivity(levpops,ne,T0,x,nH,iSpecies,lev1,lev2,emiss1)
00069 
00070       emiss = emiss1
00071 
00072    END SUBROUTINE CalcEmiss
00073 
00074    ! Gets g, e, a, and c for OI in preparation to form matrix b
00075    SUBROUTINE OI(ne,T,g,e,a,o)
00076       REAL(KIND=qPREC), INTENT(IN) :: ne, T
00077       REAL(KIND=qPREC), INTENT(OUT) :: g(5), e(5), a(5,5), o(5,5)
00078       REAL(KIND=qPREC) :: a1(5,5), erg(5,5)
00079       INTEGER :: i, j
00080 
00081       ! Statistical weights for OI
00082       g(1) = 5.
00083       g(2) = 3.
00084       g(3) = 1.
00085       g(4) = 5.
00086       g(5) = 1.
00087 
00088       ! Energy levels of OI (eV)
00089       CALL E_Values(e,erg,iOI)
00090 
00091       ! Get A-values, Baluja, K., and Zeippen, C. 1988, J. Phys B. Atom. Mol. Opt. Phys. 21, 1455.
00092       ! Subroutine A_Values gives a as a(fin,init) but we need a as
00093       ! a(init,fin) to get matrix b correct, so we take the transpose...
00094       CALL A_Values(a1,iOI)
00095       DO i=1, 5
00096          DO j=1, 5
00097             a(i,j) = a1(j,i)
00098          END DO
00099       END DO
00100 
00101       ! Collision strengths for OI, Berrington and Burke 1981, Plan. Sp. Sci. 29, 377.
00102       ! and Le Dourneuf, M. and Neabet, R. 1976, J Phys B. Atom. Molec. Phys. 9, L241.
00103       o(1,4) = (5d0/9d0)*0.266*(T**0.97)*(T**(LOG10(T**(-0.27))))
00104       o(2,4) = (3d0/9d0)*0.266*(T**0.97)*(T**(LOG10(T**(-0.27))))
00105       o(3,4) = (1d0/9d0)*0.266*(T**0.97)*(T**(LOG10(T**(-0.27))))
00106       o(1,5) = (5d0/9d0)*0.0324*T*(T**(LOG10(T**(-0.22))))
00107       o(2,5) = (3d0/9d0)*0.0324*T*(T**(LOG10(T**(-0.22))))
00108       o(3,5) = (1d0/9d0)*0.0324*T*(T**(LOG10(T**(-0.22))))
00109       o(1,2) = 0.0987*T**1.11
00110       o(2,3) = 0.0265*T**1.17
00111       o(1,3) = 0.0292*T*(T**(LOG10(T**(-0.13))))
00112       o(4,5) = 0.105*T**0.54
00113 
00114    END SUBROUTINE OI
00115 
00116    ! Gets g, e, a, and c for NII in preparation to form matrix b
00117    SUBROUTINE NII(ne,T,g,e,a,o)
00118       REAL(KIND=qPREC), INTENT(IN) :: ne, T
00119       REAL(KIND=qPREC), INTENT(OUT) :: g(5), e(5), a(5,5), o(5,5)
00120       REAL(KIND=qPREC) :: a1(5,5), erg(5,5)
00121       INTEGER :: i, j
00122 
00123       ! Statistical weights for NII
00124       g(1) = 1.
00125       g(2) = 3.
00126       g(3) = 5.
00127       g(4) = 5.
00128       g(5) = 1.
00129 
00130       ! Energy levels of NII (eV)
00131       CALL E_Values(e,erg,iNII)
00132 
00133       ! Get A-values, Bell, Hibbert, Stafford 1995, Phys. Scr. 52, 240.
00134       !               Storey and Zeippen 2000, MNRAS 312, 813.
00135       !               Froese Fischer, C. and Saha, H. 1985, Phys Scr. 32, 181.
00136       ! Subroutine A_Values gives a as a(fin,init) but we need a as
00137       ! a(init,fin) to get matrix b correct, so we take the transpose...
00138       CALL A_Values(a1,iNII)
00139       DO i=1, 5
00140          DO j=1, 5
00141             a(i,j) = a1(j,i)
00142          END DO
00143       END DO
00144 
00145       ! Collision strengths for NII, Hudson and Bell A&A 430, 725, 2005
00146       o(1,2) = 0.43
00147       o(1,3) = 0.27
00148       o(1,4) = 0.30
00149       o(1,5) = 0.035
00150       o(2,3) = 1.15
00151       o(2,4) = 0.91
00152       o(2,5) = 0.11
00153       o(3,4) = 1.51
00154       o(3,5) = 0.18
00155       o(4,5) = 0.81
00156 
00157    END SUBROUTINE NII
00158 
00159 
00160    ! Gets g, e, a, and c for SII in preparation to form matrix b
00161    SUBROUTINE SII(ne,T,g,e,a,o)
00162       REAL(KIND=qPREC), INTENT(IN) :: ne, T
00163       REAL(KIND=qPREC), INTENT(OUT) :: g(5), e(5), a(5,5), o(5,5)
00164       REAL(KIND=qPREC) :: a1(5,5), erg(5,5)
00165       INTEGER :: i, j
00166 
00167       ! Statistical weights for SII
00168       g(1) = 4.
00169       g(2) = 4.
00170       g(3) = 6.
00171       g(4) = 2.
00172       g(5) = 4.
00173 
00174       ! Energy levels of SII (eV)
00175       CALL E_Values(e,erg,iSII)
00176 
00177       ! Get A-values, Mandoza & Zeippen 1982 MNRAS 198, 111.
00178       ! Subroutine A_Values gives a as a(fin,init) but we need a as
00179       ! a(init,fin) to get matrix b correct, so we take the transpose...
00180       CALL A_Values(a1,iSII)
00181       DO i=1, 5
00182          DO j=1, 5
00183             a(i,j) = a1(j,i)
00184          END DO
00185       END DO
00186 
00187       ! Collision strengths for SII, Keenan et al 1996, MNRAS 281, 1073.
00188       o(1,2) = 2.76*(T**(-0.15))
00189       o(1,3) = 4.14*(T**(-0.15))
00190       o(1,4) = 1.17*(T**(-0.07))
00191       o(1,5) = 2.35*(T**(-0.07))
00192       o(2,3) = 7.47*(T**(-0.10))
00193       o(2,4) = 1.79*(T**(-0.17))
00194       o(2,5) = 3.00*(T**(-0.12))
00195       o(3,4) = 2.20*(T**(-0.11))
00196       o(3,5) = 4.99*(T**(-0.15))
00197       o(4,5) = 2.71*(T**(-0.14))
00198 
00199    END SUBROUTINE SII
00200 
00201 
00202    ! Forms the matrix b
00203    SUBROUTINE Loadb(b,T,ne,g,e,a,o)
00204       REAL(KIND=qPREC), INTENT(IN) :: T, ne, g(5), e(5), a(5,5), o(5,5)
00205       REAL(KIND=qPREC) :: q12,q13,q14,q15,q23,q24,q25,q34,q35,q45,q21,q31,q41,q51,q32,q42,q52,q43,q53,q54,q5,q4,q3,q2
00206       REAL(KIND=qPREC), INTENT(OUT) :: b(5,5)
00207          
00208       q12 = (8.63d-8)*o(1,2)*(EXP(-(e(2) - e(1))/(0.8614*T)))/(g(1)*(T**0.5))
00209       q13 = (8.63d-8)*o(1,3)*(EXP(-(e(3) - e(1))/(0.8614*T)))/(g(1)*(T**0.5))
00210       q14 = (8.63d-8)*o(1,4)*(EXP(-(e(4) - e(1))/(0.8614*T)))/(g(1)*(T**0.5))
00211       q15 = (8.63d-8)*o(1,5)*(EXP(-(e(5) - e(1))/(0.8614*T)))/(g(1)*(T**0.5))
00212       q23 = (8.63d-8)*o(2,3)*(EXP(-(e(3) - e(2))/(0.8614*T)))/(g(2)*(T**0.5))
00213       q24 = (8.63d-8)*o(2,4)*(EXP(-(e(4) - e(2))/(0.8614*T)))/(g(2)*(T**0.5))
00214       q25 = (8.63d-8)*o(2,5)*(EXP(-(e(5) - e(2))/(0.8614*T)))/(g(2)*(T**0.5))
00215       q34 = (8.63d-8)*o(3,4)*(EXP(-(e(4) - e(3))/(0.8614*T)))/(g(3)*(T**0.5))
00216       q35 = (8.63d-8)*o(3,5)*(EXP(-(e(5) - e(3))/(0.8614*T)))/(g(3)*(T**0.5))
00217       q45 = (8.63d-8)*o(4,5)*(EXP(-(e(5) - e(4))/(0.8614*T)))/(g(4)*(T**0.5))
00218       q21 = (8.63d-8)*o(1,2)/(g(2)*(T**0.5))
00219       q31 = (8.63d-8)*o(1,3)/(g(3)*(T**0.5))
00220       q41 = (8.63d-8)*o(1,4)/(g(4)*(T**0.5))
00221       q51 = (8.63d-8)*o(1,5)/(g(5)*(T**0.5))
00222       q32 = (8.63d-8)*o(2,3)/(g(3)*(T**0.5))
00223       q42 = (8.63d-8)*o(2,4)/(g(4)*(T**0.5))
00224       q52 = (8.63d-8)*o(2,5)/(g(5)*(T**0.5))
00225       q43 = (8.63d-8)*o(3,4)/(g(4)*(T**0.5))
00226       q53 = (8.63d-8)*o(3,5)/(g(5)*(T**0.5))
00227       q54 = (8.63d-8)*o(4,5)/(g(5)*(T**0.5))
00228       q5 = q51 + q52 + q53 + q54
00229       q4 = q41 + q42 + q43 + q45
00230       q3 = q31 + q32 + q34 + q35
00231       q2 = q21 + q23 + q24 + q25
00232 
00233       ! Matrix values before inversion (row,column)
00234       b(1,1) = ne*q15
00235       b(1,2) = ne*q25
00236       b(1,3) = ne*q35
00237       b(1,4) = ne*q45
00238       b(1,5) = -a(5,5) - ne*q5
00239       b(2,1) = ne*q14
00240       b(2,2) = ne*q24
00241       b(2,3) = ne*q34
00242       b(2,4) = -a(4,4) - ne*q4
00243       b(2,5) = a(5,4) + ne*q54
00244       b(3,1) = ne*q13
00245       b(3,2) = ne*q23
00246       b(3,3) = -a(3,3) - ne*q3
00247       b(3,4) = a(4,3) + ne*q43
00248       b(3,5) = a(5,3) + ne*q53
00249       b(4,1) = ne*q12
00250       b(4,2) = -a(2,1) - ne*q2
00251       b(4,3) = a(3,2) + ne*q32
00252       b(4,4) = a(4,2) + ne*q42
00253       b(4,5) = a(5,2) + ne*q52
00254       b(5,:) = 1.
00255 
00256    END SUBROUTINE Loadb
00257 
00258 
00259 
00260    SUBROUTINE ludcmp(b,indx)
00261       REAL(KIND=qPREC), PARAMETER :: tiny = 1d-28
00262       REAL(KIND=qPREC) :: b(5,5), indx(5), vv(5), aamax, dum, sum
00263       REAL(KIND=qPREC) :: s
00264       INTEGER :: i, j, k, imax
00265 
00266       s = 1.
00267       DO i=1, 5
00268          aamax = 0.
00269          DO j=1, 5
00270             IF(ABS(b(i,j)) > aamax) aamax = ABS(b(i,j))
00271          END DO
00272          IF(aamax == 0.) PAUSE 'singular matrix'
00273          vv(i) = 1./aamax
00274       END DO
00275       DO j=1, 5
00276          DO i=1, j-1
00277             sum = b(i,j)
00278             DO k=1, i-1
00279                sum = sum - b(i,k)*b(k,j)
00280             END DO
00281             b(i,j) = sum
00282          END DO
00283          aamax = 0.
00284          DO i=j, 5
00285             sum = b(i,j)
00286             DO k=1, j-1
00287                sum = sum-b(i,k)*b(k,j)
00288             END DO
00289             b(i,j) = sum
00290             dum = vv(i)*ABS(sum)
00291             IF(dum >= aamax) THEN
00292                imax = i
00293                aamax = dum
00294             END IF
00295          END DO
00296          IF(j /= imax) THEN
00297             DO k=1, 5
00298                dum = b(imax,k)
00299                b(imax,k) = b(j,k)
00300                b(j,k) = dum
00301             END DO
00302             s = -s
00303             vv(imax) = vv(j)
00304          END IF
00305          indx(j) = imax
00306          IF(b(j,j) == 0) b(j,j) = tiny
00307          IF(j /= 5) THEN
00308             dum = 1./b(j,j)
00309             DO i=j+1, 5
00310                b(i,j) = b(i,j)*dum
00311             END DO
00312          END IF
00313       END DO    
00314    END SUBROUTINE ludcmp
00315 
00316 
00317 
00318    SUBROUTINE lubksb(b,indx,den)
00319       REAL(KIND=qPREC) :: b(5,5), indx(5), den(5), sum
00320       INTEGER :: i, ii, j, ll
00321 
00322       ii = 0
00323       DO i=1, 5
00324          ll = indx(i)
00325          sum = den(ll)
00326          den(ll) = den(i)
00327          IF(ii /= 0) THEN
00328             DO j=ii, i-1
00329                sum = sum - b(i,j)*den(j)
00330             END DO
00331          ELSE IF(sum /= 0) THEN
00332             ii = i
00333          END IF
00334          den(i) = sum
00335       END DO
00336       DO i=5, 1, -1
00337          sum = den(i)
00338          DO j=i+1, 5
00339             sum = sum - b(i,j)*den(j)
00340          END DO
00341          den(i) = sum/b(i,i)
00342       END DO
00343    END SUBROUTINE lubksb
00344 
00345 
00346    !  Calculate the emissivity from formula (n*A*hv/(4Pi))
00347    !  Units are erg/cm^3/s/str
00348    SUBROUTINE emissivity(levpops,ne,T0,x,nH,iSpecies,lev1,lev2,emiss7)
00349       REAL(KIND=qPREC), PARAMETER :: Pi = 3.14159, k = 8.617d-5,    ! k is in eV/K
00350                      
00351                                   ! Ionization energies of neutral atoms
00352                                   EionOI = 13.62, &
00353                                   EionNI = 14.5341, &
00354                                   EionHI = 13.60, &
00355 
00356                                   ! Solar abundances from Hartigan & Morse, APJ, 660, 2007.
00357                                   AbundO = 8.82, &
00358                                   AbundN = 7.96, &
00359                                   AbundS = 7.30
00360 
00361       REAL(KIND=qPREC), INTENT(IN) :: levpops, ne, x, T0, nH
00362       INTEGER, INTENT(IN) :: lev1, lev2, iSpecies
00363       REAL(KIND=qPREC), INTENT(OUT) :: emiss7
00364       REAL(KIND=qPREC) :: kT, Abund_rat, nHII_nHI, levdens, avals(5,5), e(5), hv_ergs(5,5), 
00365                           nSpecies, nSpec, gHI, gHII,  gOI, gOII, gNI, gNII, nIon_Neut
00366 
00367       ! Statistical Weights (ground states)
00368       gHI = 2 ; gHII = 1 ; gOI = 3 ; gOII = 4 ; gNI = 4 ; gNII = 9
00369 
00370       ! kT values in eV
00371       kT = k * T0
00372 
00373       ! Calculate ratio of nHII/nHI = Xh/(1-Xh)
00374       IF(1. - x == 0.0) THEN
00375          nHII_nHI = x/(1. - x + 1d-10)
00376       ELSE
00377          nHII_nHI = x/(1. - x)
00378       END IF
00379 
00380       SELECT CASE(iSpecies)
00381          CASE(iOI)
00382             Abund_rat = AbundO - 12.
00383             nIon_Neut = nHII_nHI * (gHI*gOII/(gHII*gOI)) * EXP((EionHI - EionOI)/kT)
00384          CASE(iNII)
00385             Abund_rat = AbundN - 12.
00386             nIon_Neut = (1d0/nHII_nHI) * (gHII*gNI/(gHI*gNII)) * EXP(-(EionHI - EionNI)/kT)
00387          CASE(iSII_6716,iSII_6731)
00388             Abund_rat = AbundS - 12.
00389             nIon_Neut = 0d0           ! All S is assumed to be SII
00390       END SELECT
00391       
00392       nSpecies = nH * (10d0**Abund_rat)
00393       nSpec = nSpecies / (1d0 + nIon_Neut)
00394       levdens = levpops*nSpec
00395 
00396       CALL A_Values(avals,iSpecies)
00397       CALL E_Values(e,hv_ergs,iSpecies)
00398       emiss7 = levdens * avals(lev2,lev1) * hv_ergs(lev2,lev1) / (4d0*Pi)
00399 
00400    END SUBROUTINE emissivity
00401 
00402 
00403    ! Gives A_Values as a(final level, initial level) for given species
00404    SUBROUTINE A_Values(a, iSpecies)
00405       REAL(KIND=qPREC), INTENT(OUT) :: a(5,5)
00406       INTEGER, INTENT(IN) :: iSpecies
00407 
00408       SELECT CASE(iSpecies)
00409          CASE(iOI)
00410             a(4,5) = 1.215
00411             a(3,5) = 0d0
00412             a(2,5) = 7.60d-2
00413             a(1,5) = 2.73d-4
00414             a(3,4) = 8.93d-7
00415             a(2,4) = 1.82d-3
00416             a(1,4) = 5.63d-3
00417             a(2,3) = 1.74d-5
00418             a(1,3) = 1.20d-10
00419             a(1,2) = 8.96d-5
00420          CASE(iNII)
00421             a(4,5) = 0.923
00422             a(3,5) = 1.435d-4
00423             a(2,5) = 4.588d-2
00424             a(1,5) = 0d0
00425             a(3,4) = 3.015d-3
00426             a(2,4) = 9.819d-4
00427             a(1,4) = 1.928d-6
00428             a(2,3) = 7.46d-6
00429             a(1,3) = 1.062d-12
00430             a(1,2) = 2.08d-6
00431          CASE(iSII,iSII_6716,iSII_6731)
00432             a(4,5) = 1.03d-6
00433             a(3,5) = 0.179
00434             a(2,5) = 0.133
00435             a(1,5) = 0.225
00436             a(3,4) = 7.79d-2
00437             a(2,4) = 0.163
00438             a(1,4) = 9.06d-2
00439             a(2,3) = 3.35d-7
00440             a(1,3) = 2.60d-4
00441             a(1,2) = 8.82d-4
00442       END SELECT           
00443 
00444       a(5,5) = a(4,5) + a(3,5) + a(2,5) + a(1,5)
00445       a(4,4) = a(3,4) + a(2,4) + a(1,4)
00446       a(3,3) = a(2,3) + a(1,3)
00447       a(2:5,2) = 0d0
00448       a(4:5,3) = 0d0
00449       a(5,4) = 0d0
00450 
00451       ! Value for transitions out of level 1 = 0.0
00452       a(:,1) = 0d0
00453 
00454    END SUBROUTINE A_Values
00455   
00456 
00457    ! Gives energy levels for given species in eV, and transition energies in erg
00458    SUBROUTINE E_Values(e,erg,iSpecies)
00459       REAL(KIND=qPREC), INTENT(OUT) :: e(5), erg(5,5)
00460       INTEGER, INTENT(IN) :: iSpecies
00461       INTEGER :: i, j
00462 
00463       ! Energy levels in eV
00464       SELECT CASE(iSpecies)
00465          CASE(iOI)
00466             e(2) = 0.0196
00467             e(3) = 0.0282
00468             e(4) = 1.9689
00469             e(5) = 4.1931
00470          CASE(iNII)
00471             e(2) = 0.00604
00472             e(3) = 0.0162
00473             e(4) = 1.900
00474             e(5) = 4.056
00475          CASE(iSII,iSII_6716,iSII_6731)
00476             e(2) = 1.843
00477             e(3) = 1.847
00478             e(4) = 3.043
00479             e(5) = 3.049
00480       END SELECT
00481 
00482       e(1) = 0d0
00483       
00484       DO i=2, 5     ! initial
00485          DO j=1, 5  ! final
00486             IF (j<i) THEN
00487                erg(j,i) = (e(i) - e(j)) * (1.602d-12)
00488             ELSE
00489                erg(j,i) = 0.0
00490             END IF
00491          END DO
00492       END DO
00493 
00494       erg(:,1) = 0d0      ! Can't transition down from level 1 
00495 
00496    END SUBROUTINE E_Values
00497 
00498 
00499    ! Calculates emissivity of H-alpha from recombination and collisional
00500    ! excitation from levels 1 --> 3,4,5   
00501    SUBROUTINE alpha_emiss(ne,T0,x,nH,emissAlpha)
00502       REAL(KIND=qPREC), PARAMETER :: Halpha_hv = 3.0263d-12,  ! ergs, = 1.88889 eV, for H n= 3->2 transition
00503                                      E_thr_n3 = 12.09d0,        ! in eV, for n=1->3
00504                                      E_thr_n4 = 12.75d0,        ! in eV, for n=1->4
00505                                      E_thr_n5 = 13.056d0,       ! in eV, for n=1->5
00506                                      k_B = 8.617d-5,          ! in eV/K
00507                                      vern_a = 7.982d-11,      ! cm^3/s
00508                                      vern_b = 0.7480d0, 
00509                                      vern_T0 = 3.148d0,       ! K
00510                                      vern_T1 = 7.036d5,       ! K
00511                                      Pi = 3.14159
00512 
00513       CHARACTER(LEN = 80) :: filename
00514 
00515       INTEGER :: i,j,k                                ! counters
00516 
00517       REAL(KIND=qPREC), DIMENSION(4,25) :: coll_str   ! table of effective collision strengths (e-H)
00518  
00519       REAL(KIND=qPREC) :: recombcoef,       ! recombination rate coefficient
00520                           collexcite,       ! collisional excitation rate coefficient
00521                           emissAlpha_rec,   ! emissivity due to H-alpha recombination
00522                           emissAlpha_coll,  ! emissivity due to collisional excitation
00523                           TeLog              ! Log of Temperature (in K)
00524 
00525       REAL(KIND=qPREC), INTENT(OUT) :: emissAlpha ! Total H-alp emissivity
00526 
00527       REAL(KIND=qPREC), INTENT(IN) :: ne, T0, x, nH      ! electron dens, temp (K), ion frac, hydrogen dens
00528       REAL(KIND=qPREC) :: upsilon(3)         ! effective collision strengths
00529 
00530       !!! Recombination !!!
00531       recombcoef = vern_a / ( SQRT(T0/vern_T0) * (1d0 + SQRT(T0/vern_T0))**(1d0 - vern_b) * &
00532                    (1d0 + SQRT(T0/vern_T1))**(1d0 + vern_b) )
00533 
00534       emissAlpha_rec = Halpha_hv * ne * x*nH * recombcoef / (4d0*Pi)
00535 
00536       !!! Collisional Excitation !!!
00537       OPEN(49, file='TABLES/Anderson_Coll_Str.tab', ACTION='READ')    ! First column is T, columns 2,3,4 are
00538       READ(49, '(4(ES12.4))') coll_str                        ! log(collision strength) from level 1
00539       CLOSE(49)                                                ! to levels 3,4,5 respectively
00540 
00541       TeLog = LOG10(T0)
00542       upsilon = 0d0
00543 
00544       DO i=1, 24
00545          IF(TeLog >= coll_str(1,i) .AND. TeLog < coll_str(1,i+1)) THEN
00546 
00547             upsilon(1) = 10d0**( (TeLog - coll_str(1,i))/(coll_str(1,i+1) - coll_str(1,i)) * &
00548                          (coll_str(2,i+1) - coll_str(2,i)) + coll_str(2,i) )
00549 
00550             upsilon(2) = 10d0**( (TeLog - coll_str(1,i))/(coll_str(1,i+1) - coll_str(1,i)) * &
00551                          (coll_str(3,i+1) - coll_str(3,i)) + coll_str(3,i) )
00552 
00553             upsilon(3) = 10d0**( (TeLog - coll_str(1,i))/(coll_str(1,i+1) - coll_str(1,i)) * &
00554                          (coll_str(4,i+1) - coll_str(4,i)) + coll_str(4,i) )
00555 
00556          ELSE IF(TeLog < coll_str(1,i)) THEN       ! collision strengths are tiny at low T, so
00557             upsilon(1) = 10d0**coll_str(2,1)       ! if TeLog is below min T in the table, set
00558             upsilon(2) = 10d0**coll_str(3,1)       ! upsilon to corresponding min coll_str
00559             upsilon(3) = 10d0**coll_str(4,1)
00560          END IF
00561       END DO
00562 
00563       ! Calculate rate coefficient (cm^3/s)
00564       collexcite = 8.63d-6/(2d0*SQRT(T0)) * ( EXP(-E_thr_n3/(k_B*T0)) * upsilon(1) + &
00565                    EXP(-E_thr_n4/(k_B*T0)) * upsilon(2) + EXP(-E_thr_n5/(k_B*T0)) * upsilon(3) )
00566 
00567       emissAlpha_coll = Halpha_hv * ne * (nH - x*nH) * collexcite / (4d0*Pi)
00568 
00569       !!! Total Emission for H-alpha !!!
00570       emissAlpha = emissAlpha_rec + emissAlpha_coll
00571 
00572    END SUBROUTINE alpha_emiss
00573 
00574 END MODULE Emissions
00575 
 All Classes Files Functions Variables