Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! neqcooling.f90 is part of AstroBEAR. 00008 ! 00009 ! AstroBEAR is free software: you can redistribute it and/or modify 00010 ! it under the terms of the GNU General Public License as published by 00011 ! the Free Software Foundation, either version 3 of the License, or 00012 ! (at your option) any later version. 00013 ! 00014 ! AstroBEAR is distributed in the hope that it will be useful, 00015 ! but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 ! GNU General Public License for more details. 00018 ! 00019 ! You should have received a copy of the GNU General Public License 00020 ! along with AstroBEAR. If not, see <http://www.gnu.org/licenses/>. 00021 ! 00022 !######################################################################### 00023 00024 !============================================================================ 00025 ! BARCLAW Astrophysical Applications 00026 !============================================================================ 00027 ! Andrew Cunningham, Alexei Poludnenko, Adam Frank 00028 ! Department of Physics and Astronomy 00029 ! University of Rochester 00030 ! Rochester, NY 00031 !---------------------------------------------------------------------------- 00032 ! File: neqcooling.f90 00033 ! Type: module 00034 ! Purpose: molecular cooling and ionization (below T=10^4) 00035 ! inhomogeneous terms plus integration scheme 00036 ! Revision History: Ver. 1.0 July 2003 A. Cunningham 00037 ! Ver. 2.0 July 2003 A. Cunningham 00038 ! Added H2 cooling for high HI, density, critical density ~ 800/cc 00039 ! ----------------------------------------------------------------------- 00040 ! Features: 00041 ! H ionization and recombination for T < 10^5 (may be used to take over for bbc below 10^4 ...) 00042 ! H2 cooling for T < 10^6 00043 ! H2 dissociation recombination T < 10^5 00044 ! OI cooling using H2 as tracer for OI T < 1000 00045 ! dust cooling 00046 ! He+ ion recombination for T < 10^5 (may be used to take over for bbc below 10^4 ...) 00047 ! 00048 MODULE NEQCoolingSrc 00049 ! USE GlobalDeclarations,ONLY:cellinfo,& 00050 ! ISINFNAN,srcPrecision,& !real array (-1:32) 00051 ! muH,muH2,muHe,mue,Xmu,boltzmann,ev,amu,& !real constants 00052 ! BindH2,IonH,IonHe,gammaH,gammaH2,&!real 00053 ! rScale,nScale,TempScale,pScale,lScale,VelScale,&!real 00054 ! Lumin,mCentral,alpha,MinTemp,RunTimesc,ViscCD,ScaleGrav,ScaleCool,EOSConstants,&!real 00055 ! verbosity,& !integer array (3) 00056 ! NrHydroVars,nDim,& !integer 00057 ! iH2,iH,iHII,IH,iHe,iHeII,iHeIII,iE,iEntropy,nSpecies,nSpeciesLO,nSpeciesHI,iSpeedHI,&!interger indecies 00058 ! lH,lH2,lHII,lHe,lHeII,lHeIII,NrCons, iCooling, !logical 00059 ! USE NodeInfoDef 00060 USE DataDeclarations 00061 USE PhysicsDeclarations 00062 USE GlobalDeclarations 00063 USE SpecFun, ONLY: eone 00064 IMPLICIT NONE 00065 PRIVATE 00066 TYPE CELLINFO 00067 ! t= time at hydro step 00068 ! dt=hydro timestep 00069 REAL(KIND=qprec) :: t,dt,ScaleGravR3,precision 00070 REAL(KIND=xprec),DIMENSION(3) :: x,dx 00071 INTEGER :: iteration 00072 REAL(KIND=qprec), POINTER, DIMENSION(:) :: qScale 00073 00074 END TYPE CELLINFO 00075 INTEGER :: Verbosity(3)=0 00076 REAL(KIND=qprec), PUBLIC :: nmine=0.01, ! minimum allowed electron density 00077 Tdust=20, ! dust temp 00078 Tfloor=200, 00079 OIFrac=0.000851, ! nOatom/nHneucei 00080 metalFrac=0.01, ! nMetals/nHneucei 00081 gammac=1.666667, 00082 metal_loss ! cooling due to metal excitation 00083 00084 LOGICAL,PUBLIC :: lNeqCooling=.True., 00085 lBBC=.False., ! use BBC? 00086 lProtectNeq=.True., 00087 lchemion=.True. 00088 00089 INTEGER, PUBLIC :: NrSpecies,nSpecieslo,nSpecieshi,iH2,iH,iHII,iHe,iHeII,iHeIII 00090 LOGICAL, PUBLIC :: lH=.True.,lH2=.False.,lHII=.True.,lHe=.True.,lHeII=.False., lHeIII=.False. 00091 00092 00093 ! Delgarno McCray cooling table 00094 REAL (KIND=qPrec), PUBLIC, DIMENSION(:), ALLOCATABLE :: DMCoolingTab,HexCoolingTab,EqIonTab 00095 REAL (KIND=qPrec), DIMENSION(:,:), ALLOCATABLE :: BBC_Table,BBC_axis_vals 00096 REAL (KIND=qPrec), DIMENSION(:), ALLOCATABLE :: metal_abun,metal_states 00097 00098 REAL(KIND=qprec), PUBLIC, PARAMETER :: small = 1.d-6 00099 00100 ! Tables 00101 ! tables for analytic rate functions that can be split into temperature dependent part only 00102 ! table look ups are faster than raising to non-integer powers and divisions 00103 ! min_TAB and max_TAB are in units of log10(temp K) 00104 INTEGER,PARAMETER :: n_TAB=131072,min_TAB=0,max_TAB=7 00105 REAL(KIND=qprec),PARAMETER :: rn_TAB=n_TAB, delta_TAB=(max_TAB-min_TAB)/rn_TAB 00106 REAL(KIND=qprec), DIMENSION(n_TAB) :: H2_dust_Recomb_TAB, 00107 H2_critical_cool_TAB, 00108 dust_cool_TAB, 00109 H_recomb_TAB, 00110 H_ioniz_TAB, 00111 He_ioniz_TAB, 00112 He_recomb_TAB, 00113 HeII_ioniz_TAB, 00114 HeIII_recomb_TAB 00115 00116 00117 ! PUBLIC :: OI_cool,H2_cool,H2_critical_cool,dust_cool,H2_diss,H2_dust_Recomb,H_recomb,H_ioniz,He_ioniz,He_recomb,& 00118 ! evolve,get_mu,get_gamma,InitNeqCool,DMCoolingRate,IICoolingRate,AnalyticCoolingRate,Cool_Derivatives,Cool_Jacobian,EOS_vars,InitBBCCool,BBC_Cooling_rate,zbar_calc,InitHexCool,HexCoolingRate,EqIonVal, Addneqtracers 00119 00120 PUBLIC :: InitDMCool, DMCoolingRate, InitNeqCool, Addneqtracers, GetNEQvars, Cool_Derivatives 00121 00122 CONTAINS 00123 00124 ! Parameter LBBC tells if this module will be used with BBC ionization / cooling 00125 SUBROUTINE InitNeqCool 00126 INTEGER :: i 00127 00128 IF (LBBC) THEN 00129 CALL InitBBCCool 00130 END IF 00131 CALL InitDMCool 00132 CALL InitHexCool 00133 ! CALL EqIonInit 00134 00135 FORALL(i=1:n_TAB) H2_dust_Recomb_TAB(i) = H2_dust_Recomb(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB),Tdust) 00136 FORALL(i=1:n_TAB) H_recomb_TAB(i) = H_recomb(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)) 00137 DO i=1,n_TAB; H_ioniz_TAB(i) = H_ioniz(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)); END DO 00138 DO i=1,n_TAB; He_ioniz_TAB(i) = He_ioniz(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)); END DO 00139 DO i=1,n_TAB; HeII_ioniz_TAB(i) = HeII_ioniz(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)); END DO 00140 DO i=1,n_TAB; HeIII_recomb_TAB(i) = He_ioniz(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)); END DO 00141 FORALL(i=1:n_TAB) He_recomb_TAB(i) = He_recomb(1.d0,1.d0,10**((i-1)*delta_TAB+min_TAB)) 00142 FORALL(i=1:n_TAB) H2_critical_cool_TAB(i) = H2_critical_cool(1.d0,10**((i-1)*delta_TAB+min_TAB)) 00143 FORALL(i=1:n_TAB) dust_cool_TAB(i) = dust_cool(1.d0,10**((i-1)*delta_TAB+min_TAB)) 00144 00145 ! H_ioniz_TAB = Log10(H_ioniz_TAB) 00146 ! H_recomb_TAB = Log10(H_recomb_TAB) 00147 ! He_ioniz_TAB = Log10(He_ioniz_TAB) 00148 ! He_recomb_TAB = Log10(He_recomb_TAB) 00149 ! HeIII_recomb_TAB = Log10(HeIII_recomb_TAB) 00150 ! HeII_ioniz_TAB = Log10(HeII_ioniz_TAB) 00151 00152 END SUBROUTINE InitNeqCool 00153 00154 SUBROUTINE Addneqtracers 00155 00156 nSpecieslo=GetNextTracerIndex() !MAX(iEntropy,iE)+1 00157 00158 iH2=0.;iH=0.;iHII=0.;iHe=0.;iHeII=0.;iHeIII=0. 00159 IF(lH) THEN 00160 CALL AddTracer(iH,'H') 00161 PRINT*,'iH=',iH 00162 ELSE 00163 iH = nSpecieslo - 1 00164 END IF 00165 IF(lH2) THEN 00166 CALL AddTracer(iH2, 'H2') 00167 PRINT*,'iH2=',iH2 00168 END IF 00169 IF(lHII) THEN 00170 CALL AddTracer(iHII, 'HII') 00171 PRINT*,'iHII=',iHII 00172 END IF 00173 IF(lHe) THEN 00174 CALL AddTracer(iHe, 'He') 00175 PRINT*,'iHe=',iHe 00176 END IF 00177 IF(lHeII) THEN 00178 CALL AddTracer(iHeII, 'HeII') 00179 PRINT*,'iHeII=',iHeII 00180 END IF 00181 IF(lHeIII) THEN 00182 CALL AddTracer(iHeIII, 'HeIII') 00183 PRINT*,'iHeIII=',iHeIII 00184 END IF 00185 nSpecieshi = GetNextTracerIndex() - 1 00186 nrSpecies = nSpecieshi-nSpecieslo+1 00187 00188 END SUBROUTINE Addneqtracers 00189 00190 ! returns number density vector for each species 00191 ! if present, returns temperature in K, everything else in computational units 00192 ! NOTE : -emicro is the thermal + magnetic energy density 00193 ! -eth is thermal energy only 00194 ! -gamma = 1+P/(rho e) 00195 ! The gammas are not the same for the case of a Real Gas (Colella & Glaz 84). 00196 ! I have invoked a Real Gas EOS to provide a component of magnetic pressure to 00197 ! support the collapse of strongly cooled shocks in pure hydro sumulations 00198 SUBROUTINE EOS_vars(q,nvec,T) 00199 IMPLICIT NONE 00200 REAL(KIND=qprec), DIMENSION(1:NrHydroVars), INTENT(IN) :: q 00201 REAL(KIND=qprec), DIMENSION(1:NrHydroVars) :: qin 00202 REAL(KIND=qprec), DIMENSION(iH:nSpeciesHi), INTENT(OUT), OPTIONAL :: nvec 00203 REAL(KIND=qprec), INTENT(OUT), OPTIONAL :: T 00204 ! Internal variables 00205 REAL(KIND=qPrec) :: kein,ethin,emicroin,mein 00206 REAL(KIND=qprec), DIMENSION(iH:nSpeciesHi) :: nvecin 00207 REAL(KIND=qprec) :: muin 00208 INTEGER :: ndim_in 00209 00210 mein = 0.d0 00211 nvecin = 0.d0 00212 qin = q 00213 IF(lH2) nvecin(iH2) = qin(iH2)/muH2 00214 IF(lHII) nvecin(iHII) = qin(iHII) 00215 IF(lHe) nvecin(iHe) = qin(iHe)/muHe 00216 IF(lHeII) nvecin(iHeII) = qin(iHeII)/muHe 00217 IF(lHeIII) nvecin(iHeIII) = qin(iHeIII)/muHe 00218 IF(lH) THEN 00219 nvecin(iH) = qin(iH) 00220 ELSE 00221 ! if we are not tracking H explicitly, then it must be what is left over 00222 nvecin(iH) = (qin(1)-SUM(qin(nSpecieslo:nspecieshi)))/muH 00223 END IF 00224 nvecin = nvecin*nScale 00225 nvecin = MAX(nvecin,zero) 00226 IF(PRESENT(nvec)) nvec = nvecin(IH:nSpeciesHi) 00227 00228 ! IF(PRESENT(T) .OR. PRESENT(mu) .OR. PRESENT(eth) & 00229 ! .OR. PRESENT(emicro) .OR. PRESENT(ke)) THEN 00230 ! muin = get_mu(nvecin(IH:nSpeciesHi)) 00231 ! gammain = get_gamma(nvecin(IH:nSpeciesHi)) 00232 ! IF(PRESENT(gamma)) gamma=gammain 00233 ! IF(PRESENT(mu)) mu=muin 00234 IF(PRESENT(T)) THEN 00235 ! cooling via EOS 00236 ! kein = half*DOT_PRODUCT(qin(m_low:m_high),qin(m_low:m_high))/(qin(1)) 00237 kein = half*qin(2)**2d0/qin(1) 00238 ! IF(PRESENT(ke)) ke = kein 00239 IF(lMHD) THEN 00240 mein = half*(qin(iBx)**2d0 + qin(iBy)**2d0 + qin(iBz)**2d0) 00241 END IF 00242 ethin = qin(iE) - kein - mein 00243 00244 ! emicroin = ethin 00245 ! IF(PRESENT(emicro)) emicro = ethin 00246 ! IF(PRESENT(eth)) eth = ethin 00247 IF(PRESENT(T)) THEN 00248 T = TempScale*(gamma-1d0)*(ethin)/qin(1) 00249 END IF 00250 ! END IF 00251 END IF 00252 END SUBROUTINE EOS_vars 00253 00254 PURE FUNCTION get_gamma(nvec) 00255 REAL(KIND=qprec),INTENT(IN), DIMENSION(IH:nSpeciesHi) :: nvec 00256 REAL(KIND=qprec) get_gamma,n,ne,nH2 00257 00258 CALL nparticle(nvec,npart=n) 00259 IF(lH2) THEN 00260 nH2 = nvec(iH2) 00261 ELSE 00262 nH2=0.d0 00263 END IF 00264 IF(nH2 < small) THEN 00265 get_gamma = gammaH 00266 ELSE 00267 get_gamma = ((n-nH2)/n*gammaH/(gammaH-1.d0) + nH2/n*gammaH2/(gammaH2-1.d0)) / & 00268 ((n-nH2)/n*1.d0/(gammaH-1.d0) + nH2/n*1.d0/(gammaH2-1.d0)) 00269 END IF 00270 END FUNCTION get_gamma 00271 00272 SUBROUTINE GetNEQvars(qin,muout,nvecout) 00273 REAL(KIND=qprec),INTENT(IN) :: qin(:) 00274 REAL(KIND=qprec), INTENT(OUT) :: muout 00275 REAL(KIND=qPREC), DIMENSION(0:nSpeciesHi), INTENT(OUT) :: nvecout 00276 REAL(KIND=qprec) :: nneuc, ne, npart 00277 00278 nvecout = 0d0 00279 IF(lH2) nvecout(iH2) = qin(iH2)/muH2 00280 IF(lHII) nvecout(iHII) = qin(iHII) 00281 IF(lHe) nvecout(iHe) = qin(iHe) 00282 IF(lHeII) nvecout(iHeII) = qin(iHeII) 00283 IF(lHeIII) nvecout(iHeIII) = qin(iHeIII)/muHe 00284 IF(lH) THEN 00285 nvecout(iH) = qin(iH)/muH 00286 ELSE 00287 ! if we are not tracking H explicitly, then it must be what is left over 00288 nvecout(iH) = (qin(1)-SUM(qin(nSpecieslo:nspecieshi)))/muH 00289 END IF 00290 nvecout = nvecout*nScale 00291 nneuc = nvecout(iH2)+nvecout(iH)+nvecout(iHII)+nvecout(iHe)+nvecout(iHeII)+nvecout(iHeIII) 00292 ne = nvecout(iHII) + nvecout(iHeII) + 2d0*nvecout(iHeIII) + nmine*nneuc 00293 npart = nneuc + ne 00294 00295 ! note: we assume zero electron mass in this calculation 00296 ! ionization processes affect the mean molecular weight of the 00297 ! gas dramaticly, since energy is conserved, this reduces the 00298 ! temperature of the gas 00299 ! ie) mu ~ 0.5 for a HI gas 00300 00301 00302 muout = (nvecout(iH2)*muH2 + (nvecout(iH) + nvecout(iHII))*muH + & 00303 (nvecout(iHe) + nvecout(iHeII) + nvecout(iHeIII))*muHe)/npart 00304 00305 END SUBROUTINE GetNEQvars 00306 00307 PURE SUBROUTINE nparticle(nvec,npart,nneuc,ne) 00308 REAL(KIND=qprec), DIMENSION(IH:nSpeciesHi),INTENT(IN) :: nvec 00309 REAL(KIND=qprec), DIMENSION(0:nSpeciesHi) :: nvecin 00310 REAL(KIND=qprec), INTENT(OUT),OPTIONAL :: npart,nneuc,ne 00311 ! internal veriables 00312 REAL(KIND=qprec) :: npartin,nneucin,nein 00313 00314 nvecin=0. 00315 nvecin(IH:nSpeciesHi)=nvec 00316 00317 nneucin = nvecin(iH2)+nvecin(iH)+nvecin(iHII)+nvecin(iHe)+nvecin(iHeII)+nvecin(iHeIII) 00318 00319 ! # of electrons 00320 nein = nvecin(iHII) + nvecin(iHeII) + 2.0*nvecin(iHeIII)+nmine*nneucin 00321 ! # of particles in gas, including electrons 00322 npartin = nneucin + nein 00323 ! # of neuclei (either molecular, neutral, or ionized) 00324 nneucin = nneucin + nvecin(iH2) 00325 00326 IF(PRESENT(ne)) ne = nein 00327 IF(PRESENT(nneuc)) nneuc = nneucin 00328 IF(PRESENT(npart)) npart = npartin 00329 END SUBROUTINE nparticle 00330 00331 ! computes energy source term from H2 radiative cooling 00332 ! See : Lepp & Schull 1983 (1983ApJ...270..578L), 00333 ! O'Sullivan & Ray 2000 (2000A&A...363..355O) 00334 PURE FUNCTION H2_cool(nH2, nH, T) !erg*cm^-3*s^-1 00335 REAL(KIND=qprec),INTENT(IN)::nH, nH2, T 00336 REAL(KIND=qprec) H2_cool,Tin,heat_H2_cool 00337 00338 IF(nH2 < small) THEN 00339 H2_cool=0. 00340 ELSE 00341 ! this function is only good to 10^6 00342 Tin = min(1.d6,T) 00343 heat_H2_cool = (Lvh(MinTemp)/(1.+Lvh(MinTemp)/Lvl(nH2,nH,MinTemp)) + Lrh(MinTemp)/(1.+Lrh(MinTemp)/Lrl(nH2,nH,MinTemp))) 00344 H2_cool = nH2*(Lvh(Tin)/(1.+Lvh(Tin)/Lvl(nH2,nH,Tin)) + Lrh(Tin)/(1.+Lrh(Tin)/Lrl(nH2,nH,Tin)) - heat_H2_cool) 00345 END IF 00346 CONTAINS 00347 PURE FUNCTION Lvh(T) ! ergs/s 00348 REAL(KIND=qprec),INTENT(IN)::T 00349 REAL(KIND=qprec) Lvh 00350 ! Smith & Rosen 2003 correct error of Lepp & Schull 1.1d-18, NOT 1.1d-13 00351 Lvh = 1.10d-18*exp(-6744./T) 00352 END FUNCTION Lvh 00353 PURE FUNCTION Lvl(nH2,nh,T) ! ergs/s 00354 REAL(KIND=qprec),INTENT(IN)::nh,nH2,T 00355 REAL(KIND=qprec) Lvl 00356 Lvl = 8.18d-13*(nH*kH(T)+nH2*kH2(T)) 00357 END FUNCTION Lvl 00358 PURE FUNCTION Qfunc(nH2,nH) 00359 REAL(KIND=qprec),INTENT(IN)::nH2,nH 00360 REAL(KIND=qprec) Qfunc!,nH2in,nHin 00361 !nH2in = max(0.d0,nH2) 00362 !nHin = max(0.d0,nH) 00363 Qfunc = nH2**.77 + 1.2*nH**.77 00364 END FUNCTION Qfunc 00365 PURE FUNCTION kh2(T) ! cm^3/s 00366 REAL(KIND=qprec),INTENT(IN)::T 00367 REAL(KIND=qprec) kh2 00368 ! deltaE = v=0->1 energy gap of H2 in eV 00369 ! nH is not zero. 00370 REAL(KIND=qprec),PARAMETER :: deltaE = 0.785 00371 kh2 = 1.45d-12*sqrt(T)*exp(-4.2*deltaE*ev/(boltzmann*(T+1190.))) 00372 END FUNCTION kh2 00373 ! piecewise functions: 00374 PURE FUNCTION kh(T) ! cm^3/s 00375 REAL(KIND=qprec),INTENT(IN)::T 00376 REAL(KIND=qprec) kh 00377 IF(T > 1635.) THEN 00378 kh = 1.0d-12*sqrt(T)*exp(-1000./T) 00379 ELSE 00380 kh = 1.4e-13*exp(T/125. - (T/577.)**2) 00381 END IF 00382 END FUNCTION kh 00383 00384 PURE FUNCTION Lrh(T) ! ergs/s 00385 REAL(KIND=qprec),INTENT(IN)::T 00386 REAL(KIND=qprec) Lrh,x 00387 IF(T > 1087.) THEN 00388 Lrh = 3.90d-19*exp(-6118.d0/T) 00389 ELSE 00390 x = log10(T/1.d4) 00391 Lrh = 10**(-19.24d0 + .0474d0*x - 1.247d0*x**2) 00392 END IF 00393 END FUNCTION Lrh 00394 00395 PURE FUNCTION Lrl(nH2,nH,T) ! ergs/s 00396 REAL(KIND=qprec),INTENT(IN):: nH2, nH, T 00397 REAL(KIND=qprec) Lrl,x 00398 IF(T > 4031.) THEN 00399 Lrl = 1.38d-22*exp(-9243.d0/T) 00400 ELSE 00401 x = log10(T/1.d4) 00402 Lrl = 10**(-22.90d0 - .553d0*x - 1.148d0*x**2) 00403 END IF 00404 Lrl = LrL*Qfunc(nH2,nH) 00405 END FUNCTION Lrl 00406 END FUNCTION H2_cool 00407 00408 PURE FUNCTION H2_critical_cool(nH2, T) !erg*cm^-3*s^-1 00409 REAL(KIND=qprec) :: H2_critical_cool, beta, exp_beta 00410 REAL(KIND=qprec),INTENT(IN)::nH2, T 00411 INTEGER :: j,v 00412 REAL(KIND=qprec), PARAMETER :: hc=6.6260755d-27*2.99792458d10 00413 !Rotational wavelengths & einstien coefficients for pure rotational v=0, S branch as function of J (rotation number) 00414 REAL(KIND=qprec), PARAMETER, DIMENSION(0:15) :: 00415 A=(/2.94d-11,4.76d-10,2.76d-9,9.84d-9,2.64d-8,5.88d-8,1.14d-7,2.00d-7, 00416 3.24d-7,4.90d-7,7.03d-7,9.64d-7,1.27d-6,1.62d-6,2.00d-6,2.41d-6 /), 00417 lam=(/28.2240d-4,17.0378d-4,12.2806d-4,9.6662d-4,8.0255d-4,6.9087d-4,6.1056d-4,5.5045d-4, 00418 5.0403d-4,4.6720d-4,4.3725d-4,4.1223d-4,3.9069d-4,3.7150d-4,3.5374d-4,3.3664d-4 /) 00419 ! I only include three brightest rotation-vibration lines 00420 ! rotation-vibration lines for v=1->0 as function of J (rotation number) 00421 ! S-branch 00422 REAL(KIND=qprec), PARAMETER, DIMENSION(0:1) :: 00423 Av1_0S = (/2.53d-7,3.47d-7 /), 00424 lamv1_0S = (/2.2232d-4,2.1217d-4 /), 00425 ! Q-branch 00426 Av1_0Q = (/zero,4.29d-7 /),& 00427 lamv1_0Q = (/one,2.4066d-4 /) 00428 ! 00429 beta=hc/(boltzmann*T) 00430 00431 H2_critical_cool = 0 00432 ! add up v=0 lines 00433 v=0 00434 DO j=0,15 00435 H2_critical_cool = H2_critical_cool + hc/(4.d0*PI*lam(j))*A(j)*fv(v)*fj(j+2) 00436 END DO 00437 ! add up v=1 lines 00438 DO j=0,1 00439 v=1 00440 H2_critical_cool = H2_critical_cool + hc/(4.d0*PI*lamv1_0S(j))*Av1_0S(j)*fv(v)*fj(j+2) 00441 H2_critical_cool = H2_critical_cool + hc/(4.d0*PI*lamv1_0Q(j))*Av1_0Q(j)*fv(v)*fj(j) 00442 END DO 00443 H2_critical_cool = nH2*H2_critical_cool 00444 00445 CONTAINS 00446 ! boltzmann population faction for rotational state j 00447 ! valid fo T > ~ 1000 K 00448 PURE FUNCTION fj(j) 00449 REAL(KIND=qprec) :: fj,x 00450 INTEGER,INTENT(IN) :: j 00451 00452 x=beta*Bv(v) 00453 ! limit fj to one. if T < 1000 K or so, apoximating patition function sum 00454 ! as integal beaks down and can lead to unphysical fractions 00455 fj = MIN(gs(j)*(2*j+1)*half*x*exp(-x*j*(j+1)),one) 00456 END FUNCTION fj 00457 00458 ! boltzmann population faction for vibrational state v 00459 PURE FUNCTION fv(v) 00460 REAL(KIND=qprec) :: fv,x 00461 REAL(KIND=qprec),PARAMETER :: we=4396.5011d0 00462 INTEGER,INTENT(IN) :: v 00463 00464 x=beta*we 00465 fv = (one-exp(-x))*exp(-x*v) 00466 END FUNCTION fv 00467 00468 ! degeneracy of rotational state j 00469 PURE FUNCTION gs(j) 00470 INTEGER,INTENT(IN) :: j 00471 INTEGER :: gs 00472 ! odd -> 3, even->one 00473 gs = MAX(MOD(j,2)*3,1) 00474 END FUNCTION gs 00475 00476 PURE FUNCTION Bv(v) 00477 INTEGER,INTENT(IN) :: v 00478 REAL(KIND=qprec) Bv 00479 REAL(KIND=qprec),PARAMETER ::alphae=2.9263d0,Be=60.7864d0 00480 Bv = Be-alphae*(REAL(v,qprec)+half) 00481 END FUNCTION Bv 00482 END FUNCTION H2_critical_cool 00483 00484 PURE FUNCTION H2_critical_cool_table(nH2, T) !erg*cm^-3*s^-1 00485 REAL(KIND=qprec),INTENT(IN)::nH2, T 00486 REAL(KIND=qprec) :: H2_critical_cool_table,Tin,inc 00487 INTEGER :: lower,upper ! the lower,upper interpolation point 00488 00489 Tin = log10(T) 00490 inc = (Tin-min_TAB)/delta_TAB 00491 lower = MAX(MIN(INT(inc),n_TAB),1) 00492 inc = inc-REAL(lower,qprec) 00493 upper = MAX(MIN(lower + 1,n_TAB),1) 00494 H2_critical_cool_table = nH2*(H2_critical_cool_TAB(lower)+(H2_critical_cool_TAB(upper)-H2_critical_cool_TAB(lower))*inc) 00495 END FUNCTION H2_critical_cool_table 00496 00497 ! reutuns nH/nH_critical, where ncritical is the critical density for H2 cooling 00498 ! to reach an equillibrium level population & therefore different cooling funciton 00499 ! For citical density of H2 lines see: Mandy & Martin (1993 ApJS...86..199M) 00500 ! Note that critical density (nCritical) is defined here at temperature TCritical. 00501 PURE FUNCTION CriticalDenH2(nH,T) 00502 REAL (KIND=qPrec), INTENT(IN) :: nH,T 00503 REAL (KIND=qprec) :: CriticalDenH2,x,ncr 00504 00505 ! x = n/nCritical, density as fraction of critical density 00506 x = log10(T/1.d4) 00507 ! HI citical density from Lepp & Schull 00508 ncr = 10**(4.d0-0.416d0*x-0.327d0*x**2) 00509 ! limit critical density to one particle per cc 00510 CriticalDenH2 = nH/MAX(ncr,one) 00511 END FUNCTION CriticalDenH2 00512 00513 00514 ! cooling contribution from H atoms "sticking" to dust 00515 ! see: Smith & Rosen 2003 (2003MNRAS.339..133S) 00516 ! Hollembach & Mckee 1989 00517 PURE FUNCTION dust_cool(n,T) !erg*cm^-3*s^-1 00518 REAL(KIND=qprec),INTENT(IN)::n, T 00519 REAL(KIND=qprec) dust_cool, Tin 00520 ! NOTE: C=1.66044d-24**1.5*4.9d-33 00521 REAL(KIND=qprec), PARAMETER :: C=4.9d-33 00522 00523 IF(T < 10.d0) THEN 00524 dust_cool = 0.d0 00525 ELSE 00526 dust_cool = n**2*3.8d-33*sqrt(T)*(T-Tdust)*(one-0.8d0*exp(-75.d0/T)) 00527 END IF 00528 END FUNCTION dust_cool 00529 00530 PURE FUNCTION dust_cool_table(n,T) !erg*cm^-3*s^-1 00531 REAL(KIND=qprec),INTENT(IN)::n, T 00532 REAL(KIND=qprec) :: dust_cool_table,Tin,inc 00533 INTEGER :: lower,upper ! the lower,upper interpolation point 00534 00535 Tin = log10(T) 00536 inc = (Tin-min_TAB)/delta_TAB 00537 lower = MAX(MIN(INT(inc),n_TAB),1) 00538 inc = inc-REAL(lower,qprec) 00539 upper = MAX(MIN(lower + 1,n_TAB),1) 00540 dust_cool_table = n**2*(dust_cool_TAB(lower)+(dust_cool_TAB(upper)-dust_cool_TAB(lower))*inc) 00541 END FUNCTION dust_cool_table 00542 00543 00544 ! Delgarno & McCray cooling interoplation 00545 FUNCTION DMCoolingRate(temp) 00546 REAL(KIND=qPrec) :: temp,DMCoolingRate 00547 REAL (KIND=qPrec) :: tpos, dtpos 00548 INTEGER itpos1,itpos2 00549 00550 tpos=(LOG10(Temp)-2)*10+1.0 ! Position on the cooling curve 00551 itpos1=MIN(60,MAX(1,INT(tpos))) ! Ensure that the lower bounding value is 00552 ! between 0 and 60 00553 dtpos=min(one, max(zero, tpos-REAL(itpos1,qPrec))) ! Interpolation increment 00554 itpos2=MIN(61,itpos1+1) ! Upper bounding value 00555 00556 ! Perform the interpolation on the cooling curve and find the cooling rate 00557 DMCoolingRate = DMCoolingTab(itpos1)+(DMCoolingTab(itpos2)-DMCoolingTab(itpos1))*dtpos 00558 DMCoolingRate = 10.d0**(DMCoolingRate) 00559 END FUNCTION DMCoolingRate 00560 00561 PURE FUNCTION EqIonVal(Temp) ! erg cm^3 sec^-1 00562 ! Interface declarations 00563 REAL (KIND=qPrec) ::EqIonVal 00564 REAL (KIND=qPrec), INTENT(IN) :: Temp 00565 ! Provide the value of the cooling rate (and heating rate) at a given 00566 ! temperature 00567 REAL (KIND=qPrec) :: tpos, dtpos 00568 INTEGER itpos1,itpos2 00569 00570 tpos=(LOG10(Temp)-2.0)*10.0+1.0 ! Position on the cooling curve 00571 itpos1=MIN(60,MAX(1,INT(tpos))) ! Ensure that the lower bounding value is 00572 ! between 1 and 60 00573 dtpos=tpos-REAL(itpos1,qPrec) ! Interpolation increment 00574 itpos2=MIN(61,itpos1+1) ! Upper bounding value 00575 00576 ! Perform the interpolation on the cooling curve and find the cooling rate 00577 IF(Temp .GT. 1.d4) THEN 00578 EqIonVal = (EqIonTab(itpos1)+(EqIonTab(itpos2)-EqIonTab(itpos1))*dtpos) 00579 ELSE 00580 EqIonVal = 5.0d3 00581 END IF 00582 END FUNCTION EqIonVal 00583 00584 FUNCTION BBC_Cooling_rate(X,Tin) 00585 00586 REAL(Kind=qPrec), PARAMETER :: Xmin=.01,Tmin=4.0,Xmax=.9,Tmax=8.0 00587 REAL(Kind=qPrec) :: BBC_Cooling_rate 00588 REAL(Kind=qPrec), Intent(IN) :: X,Tin !Ionization fraction and Temperature 00589 ! X = ne/nhtot 00590 REAL(Kind=qPreC) :: inc1,inc2,inc3,slope1,slope2,slope3 00591 REAL(Kind=qPrec) :: inter1,inter2,inter3,zb,T 00592 INTEGER :: Tlower,Tupper,Zlower,Zupper 00593 INTEGER ,Dimension(2):: nInputs 00594 INTEGER :: I,J 00595 00596 00597 CALL zbar_calc(T,X,zb)!,metal_abun,metal_states) 00598 T = log10(Tin) 00599 IF ((T .LE. Tmin) .OR. (T .GE. Tmax)) THEN 00600 BBC_Cooling_rate = 0.0 00601 RETURN 00602 END IF 00603 00604 !!!!!Find upper and lower bound for X and T for interpolating 00605 nInputs = shape(BBC_Table) 00606 00607 DO I=1,nInputs(2) 00608 IF (zb .GE. BBC_axis_vals(2,I)) THEN 00609 Zlower = I 00610 ELSE 00611 EXIT 00612 END IF 00613 END DO 00614 DO I=1,nInputs(1) 00615 If (T .GE. BBC_axis_vals(1,I)) THEN 00616 Tlower = I 00617 ELSE 00618 EXIT 00619 END IF 00620 END DO 00621 Zupper = Zlower + 1 00622 Tupper = Tlower + 1 00623 !!!!!!!!!!!!! 00624 00625 !!!!!!!!! 3 seperate interpolations for 2-D grid BBC Table 00626 inc1 = (zb - BBC_axis_vals(2,Zlower)) / (BBC_axis_vals(2,ZUpper)-BBC_axis_vals(2,Zlower)) 00627 slope1 = BBC_Table(Tupper,Zupper) - BBC_Table(Tupper,Zlower) 00628 inter1 = BBC_Table(Tupper,Zlower) + slope1*inc1 00629 00630 inc2 = (zb - BBC_axis_vals(2,Zlower)) / (BBC_axis_vals(2,ZUpper)-BBC_axis_vals(2,Zlower)) 00631 slope2 = BBC_Table(Tlower,Zupper) - BBC_Table(Tlower,Zlower) 00632 inter2 = BBC_Table(Tlower,Zlower) + slope2*inc2 00633 00634 inc3 = T - BBC_axis_vals(1,Tlower) 00635 slope3 = (inter1 - inter2) / (BBC_axis_vals(1,Tupper)-BBC_axis_vals(1,Tlower)) 00636 inter3 = inter2 + slope3*inc3 00637 BBC_Cooling_rate = inter3 00638 !!!!!!!!!!!!!!!!!!!!!! 00639 00640 END FUNCTION BBC_Cooling_rate 00641 00642 SUBROUTINE InitBBCCool 00643 REAL, PARAMETER :: qPrec = Selected_Real_Kind(14,32) 00644 INTEGER :: nInput1,nInput2 00645 INTEGER :: iErr,I,num_metals 00646 00647 OPEN(86,IOSTAT=iErr,file='TABLES/BBC.tab') 00648 READ(86,*) nInput1,nInput2 00649 ALLOCATE(BBC_Table(nInput1,nInput2)) 00650 ALLOCATE(BBC_axis_vals(2,MAX(nInput1,nInput2))) 00651 BBC_axis_vals =0.0 00652 READ(86,*) BBC_axis_vals(1,1:nInput1) 00653 READ(86,*) BBC_axis_vals(2,1:nInput2) 00654 READ(86,*) BBC_Table 00655 CLOSE(86) 00656 BBC_Table = BBC_Table * 1.0e-23 00657 OPEN(111,file='./BBC_metal.abun') 00658 READ(111,*) num_metals 00659 ALLOCATE(metal_abun(num_metals)) 00660 ALLOCATE(metal_states(num_metals)) 00661 DO I=1,num_metals 00662 READ(111,*) metal_abun(I),metal_states(I) 00663 END DO 00664 metal_abun = 10**(metal_abun-12) !scales to units where pure H has an abun of 12 00665 CLOSE(111) 00666 00667 00668 END SUBROUTINE InitBBCCool 00669 00672 SUBROUTINE InitDMCool 00673 INTEGER :: i,nInput 00674 REAL :: dummy 00675 00676 IF(ALLOCATED(DMCoolingTab)) RETURN 00677 OPEN(UNIT=85,FILE='TABLES/DMcooling.tab') 00678 READ(85,*) nInput 00679 ALLOCATE(DMCoolingTab(nInput)) 00680 DO i=1,nInput 00681 READ(85,*) dummy, DMCoolingTab(i) 00682 END DO 00683 CLOSE(85) 00684 END SUBROUTINE InitDMCool 00685 00686 SUBROUTINE EqIonInit 00687 ! Initialize the cooling routine and load the cooling table from a file 00688 REAL (KIND=qPrec) :: Temp 00689 INTEGER i,nInput,iErr 00690 ! 00691 OPEN(UNIT=85,IOSTAT=iErr,ERR=10,FILE='TABLES/EqIon.tab', & 00692 STATUS='old',FORM='formatted') 00693 10 IF (iErr /= 0) THEN 00694 PRINT *,'!!! src FATAL: Error opening file "TABLES/EqIon.tab"' 00695 STOP 00696 END IF 00697 ! 00698 READ(85,*) nInput 00699 ALLOCATE(EqIonTab(nInput),STAT=iErr) 00700 IF (iErr/=0) THEN 00701 PRINT *,'!!! src FATAL: Cannot allocate EqIon array' 00702 STOP 00703 END IF 00704 DO i=1,nInput 00705 READ(UNIT=85,FMT=*,IOSTAT=iErr,END=20) Temp,EqIonTab(i) 00706 20 IF (iErr /= 0) THEN 00707 PRINT *,'!!! src FATAL: EOF encountered while reading file ', & 00708 '"TABLES/EqIon.tab"' 00709 STOP 00710 END IF 00711 END DO 00712 CLOSE(85) 00713 ! Convert cooling from log to linear 00714 ! EqIonTab(1:nInput)=10.0d0**EqIonTab(1:nInput) 00715 END SUBROUTINE EqIonInit 00716 00717 00718 !!!! Cooling due to non-equilibrium hydrogen excitation 00719 SUBROUTINE InitHexCool 00720 REAL(KIND=qPREC) :: dummy 00721 INTEGER :: i,nInput 00722 00723 IF(ALLOCATED(HexCoolingTab)) RETURN 00724 OPEN(UNIT=85, FILE='TABLES/Hexcooling.tab') 00725 READ(85,*) nInput 00726 ALLOCATE(HexCoolingTab(nInput)) 00727 DO i=1, nInput 00728 READ(85,*) dummy, HexCoolingTab(i) 00729 END DO 00730 CLOSE(85) 00731 END SUBROUTINE InitHexCool 00732 00733 PURE FUNCTION HexCoolingRate(Temp) ! erg cm^3 sec^-1 00734 ! Interface declarations 00735 REAL (KIND=qPrec) :: HexCoolingRate 00736 REAL (KIND=qPrec), INTENT(IN) :: Temp 00737 ! Provide the value of the cooling rate (and heating rate) at a given 00738 ! temperature 00739 REAL (KIND=qPrec) :: tpos, dtpos 00740 INTEGER itpos1,itpos2 00741 00742 tpos=(LOG10(Temp)-2.0)*10.0+1.0 ! Position on the cooling curve 00743 itpos1=MIN(60,MAX(1,INT(tpos))) ! Ensure that the lower bounding value is 00744 ! between 1 and 60 00745 dtpos=tpos-REAL(itpos1,qPrec) ! Interpolation increment 00746 itpos2=MIN(61,itpos1+1) ! Upper bounding value 00747 00748 ! Perform the interpolation on the cooling curve and find the cooling rate 00749 IF(Temp .GT. 2.e2) THEN 00750 HexCoolingRate = (HexCoolingTab(itpos1)+(HexCoolingTab(itpos2)-HexCoolingTab(itpos1))*dtpos) 00751 HexCoolingRate = 10**HexCoolingRate 00752 ELSE 00753 HexCoolingRate = 0. 00754 END IF 00755 END FUNCTION HexCoolingRate 00756 00757 00758 ! Oxygen cooling from OI excitation due to O-H collisions 00759 ! see: Launay & Roueff (1977A&A....56..289L) 00760 ! 00761 ! This includes 63 and 44 micrometer lines due to collision with hydrogens. The 00762 ! dominant term is from the 63 micrometer line. 00763 ! 00764 ! OI electron collisions are only important for electron fraction > ~1e-4 00765 ! Federman & Shipsey (ApJ...269..791F). At temperatures high enough to attain 00766 ! significant electron fraction, H2 will have dissociated, so i think we can neglect 00767 ! OI interaction w/ e-. 00768 ! 00769 ! Federman & Shipsey (ApJ...269..791F) also say that H2 - OI collisions results 00770 ! in no line cooling but rater an OH formation. 00771 ! 00772 PURE FUNCTION OI_cool(nH2,nH,T) 00773 REAL(KIND=qprec),INTENT(IN)::nH2,nH,T 00774 REAL(KIND=qprec) OI_cool,Ttab 00775 REAL(KIND=qprec), DIMENSION(0:10), PARAMETER :: 00776 CoolTab = (/0.d0,1.91d0,10.7d0,21.3d0,31.5d0,41.1d0,50.0d0,58.2d0,65.8d0,72.9d0,79.7d0/)*1.d-25 00777 REAL(KIND=qprec), PARAMETER :: deltaT = 100, 00778 minT = 100., 00779 maxT = 1000. 00780 INTEGER(KIND=qprec) :: lo,hi,Tlo,Thi 00781 00782 IF(T <= minT) THEN 00783 OI_cool=0. 00784 RETURN 00785 ELSE IF(T < maxT) THEN 00786 Ttab = min(T,maxT) 00787 lo = INT((T-minT)/deltaT) 00788 Tlo = minT + deltaT*lo 00789 hi = lo+1 00790 Thi = minT + deltaT*hi 00791 OI_cool = CoolTab(lo) + (CoolTab(hi)-CoolTab(lo))/(Thi-Tlo)*(T-Tlo) 00792 ! if we are above 1000k, extrapolate out as sqrt(T) 00793 ELSE 00794 OI_cool = CoolTab(10)*sqrt(T/maxT) 00795 END IF 00796 ! number of O atom = OIFrac*nHneucei 00797 ! Ionization potential of OI is about the same as that of H. 00798 ! fraction of H atoms (either atomic or molecular) traces OI abundance 00799 ! OI_cool = (OIFrac*nneuclei)*(nH+2.*nH2)/nHneuclei*nH*OI_cool 00800 OI_cool = OIFrac*(nH+2.*nH2)*nH*OI_cool 00801 END FUNCTION OI_cool 00802 00803 00804 ! computes H2 dissociation / resociation source term 00805 ! recombination occurs only on dust grain surfaces 00806 ! This considers collisions with H, HE, H2, and e-, no H+ collision 00807 ! is considered. H+ fractions (ane e- for that matter) are important 00808 ! only after all of the H2 has been dissociated anyway, I think .... 00809 ! See : O'Sullivan & Ray 2000 (2000A&A...363..355O) (implementation info) 00810 ! Lim et al. (2002 RAS,MNRAS 335,817-824) (another implementation) 00811 ! + references listed in sub-functions 00812 PURE FUNCTION H2_diss(nH2, nH, nHe, ne, T) !cm^-3*s^-1 00813 REAL(KIND=qprec),INTENT(IN)::nH2, nH, nHe, ne, T 00814 REAL(KIND=qprec) H2_diss,n,Tin 00815 ! collisional partner reference indecies 00816 INTEGER H,H2,He,e 00817 PARAMETER(H=1,H2=2,He=3,e=4) 00818 H2_diss=zero 00819 IF(T<50.) RETURN 00820 n = nH2 + nH + nHe + ne 00821 ! this function can cause nasty divergence at high temperatures,limit to 10k 00822 Tin = min(1.d5,T) 00823 ! NOTE:: low temp (T<~300) rates can cause divide by zero in kd function 00824 IF(nH2 <= small .OR. Tin < 1000.d0) THEN 00825 H2_diss=zero 00826 ELSE 00827 H2_diss = nH2*kd(H2,Tin) 00828 IF(nH>small) H2_diss = H2_diss + nH*kd(H,Tin) 00829 IF(nHe>small) H2_diss = H2_diss + nHe*kd(He,Tin) 00830 IF(ne>small) H2_diss = H2_diss + ne*kd(e,Tin) 00831 H2_diss = nH2*H2_diss 00832 END IF 00833 CONTAINS 00834 00835 ! interpolates between high and low density limites of 00836 ! collisional dissociation rate for H2 + partner collision 00837 ! i = index for a particular collision partner 00838 ! H = 1 00839 ! H2 = 2 00840 ! if both limits are avalible -- otherwise it uses low density 00841 ! limit (an underestimate of the dissociation rate 00842 ! 00843 ! For high density H2 & H collisions I use the rates of Lepp & Schull 00844 ! For low densuty H collisions I use the rate of Dove & Mandy 00845 ! For low density H2 collisions I use the rate of Martin & Keogh 00846 ! For e- and He collisions I use the low density limit only (no high linit avalible) 00847 ! from Martin & Keogh. These are lower limits on the actual rate. 00848 ! The interpolation scheme between the limits is from Lepp & Schull 00849 ! 00850 ! see: Lim et al. (2002 RAS,MNRAS 335,817-824), Appendix 00851 PURE FUNCTION kd(i,T) 00852 REAL(KIND=qprec),INTENT(IN)::T 00853 INTEGER,INTENT(IN):: i 00854 REAL(KIND=qprec):: kd,kdL,kdH,ncr,x 00855 00856 kdL=0;kdH=0;kd=0;ncr=0 00857 00858 SELECT CASE(i) 00859 CASE(H) 00860 kdH = kdH_LS(i,T) 00861 kdL = kdL_DM(i,T) 00862 !kdL = kdL_LS(i,T) 00863 x = log10(T/1.d4) 00864 ncr = 10**(4.d0-0.416d0*x-0.327d0*x**2) 00865 CASE(H2) 00866 kdH = kdH_LS(i,T) 00867 kdL = kdL_MK(i,T) 00868 !kdL = kdL_LS(i,T) 00869 x = log10(T/1.d4) 00870 ncr = 10**(4.13d0-0.986d0*x+0.119d0*x**2) 00871 CASE(He,e) 00872 kd = kdL_MK(i,T) 00873 RETURN 00874 END SELECT 00875 ! interpolation function 00876 kd = 10**(log10(kdH) - log10(kdH/kdL)/(one+n/MAX(ncr,small))) 00877 END FUNCTION kd 00878 00879 ! collisional dissociation rate for H2 + partner collision 00880 ! i = index for a particular collision partner 00881 ! H = 1 00882 ! H2 = 2 00883 ! see: Lepp & Schull (1983 ApJ) 00884 ! NOTE: these parameters are only good in high density limit n >~ 1e6 00885 ! for lower densities, these are upper bounds on the dissociation rate 00886 PURE FUNCTION kdH_LS(i,T) 00887 REAL(KIND=qprec),INTENT(IN)::T 00888 INTEGER,INTENT(IN):: i 00889 REAL(KIND=qprec)::kdH_LS 00890 00891 kdH_LS = 0 00892 SELECT CASE(i) 00893 CASE(H) 00894 kdH_LS = 3.52e-9*exp(-4.39d4/T) 00895 CASE(H2) 00896 kdH_LS = 5.48e-9*exp(-5.30d4/T) 00897 END SELECT 00898 END FUNCTION kdH_LS 00899 00900 ! collisional dissociation rate for H2 + partner collision 00901 ! i = index for a particular collision partner 00902 ! H = 1 00903 ! H2 = 2 00904 ! He = 3 00905 ! e- = 4 00906 ! see: Martin & Keogh 1998 (1998ApJ...499..793M) 00907 ! NOTE: these parameters are only good in low density limit n <~ 1e3 00908 ! for larger densities, these are lower bounds on the dissociation rate 00909 ! 00910 ! For low density H2-H2 collisions these are more accurate than those given by Lepp & Schull 00911 ! according to Lim et al. (2002 RAS,MNRAS 335,817-824), Appendix 00912 PURE FUNCTION kdL_MK(i,T) 00913 REAL(KIND=qprec),INTENT(IN)::T 00914 INTEGER,INTENT(IN):: i 00915 REAL(KIND=qprec)::kdL_MK,Tin 00916 REAL(KIND=qprec), DIMENSION(4) :: a,b,c,Eo,g,mu 00917 ! d = constant to give units in cm^3 S^-1 00918 ! k = boltzman constant in Martin & Keogh's screwed up unit scale 00919 REAL(KIND=qprec) d, k 00920 ! NOTE: g = Gamma(b+1) 00921 PARAMETER(g = (/3.60358,71.1182,1.70622,1.04386/)) 00922 PARAMETER(a = (/54.1263,40.1008,4.8152,11.2474/)) 00923 PARAMETER(b = (/2.5726,4.6881,1.8208,1.0948/)) 00924 PARAMETER(c = (/3.45,2.1347,-0.9459,2.3182/)) 00925 PARAMETER(Eo = (/0.168,0.1731,0.4146,0.3237/)) 00926 PARAMETER(mu = (/muH*muH2/(muH+muH2),muH2*muH2/(muH2+muH2),muHe*muH2/(muHe+muH2),mue*muH2/(mue+muH2)/)) 00927 PARAMETER(d = 1.849d-22) 00928 PARAMETER(k = 3.167d-6) 00929 ! this function is only good to 10^5 according to Martin & Keogh 00930 Tin = min(1.d5,T) 00931 ! collisions with He are very strong above 10^5 and this number will overflow, limit to ~10^5 00932 if(i == 3) then 00933 Tin = min(1.2d5,T) 00934 else 00935 ! for species other than He, I find this function is similar to others 00936 ! below 10^6 00937 Tin = min(1.d6,T) 00938 end if 00939 kdL_MK = d*sqrt(8.d0*k*Tin/(PI*mu(i)*AMU))*a(i)*(k*Tin)**(b(i)-1.d0)*g(i)*exp(-Eo(i)/(k*Tin))/((1.+c(i)*k*Tin)**(b(i)+1.d0)) 00940 END FUNCTION kdL_MK 00941 00942 ! collisional dissociation rate for H2 + partner collision 00943 ! i = index for a particular collision partner 00944 ! H = 1 00945 ! see: Dove & Mandy (1986) 00946 ! NOTE: these parameters are only good in low density limit n <~ 1e3 00947 ! for larger densities, these are lower bounds on the dissociation rate 00948 ! For low density H-H2 collisions these are more accurate than those given by Lepp & Schull 00949 ! according to Lim et al. (2002 RAS,MNRAS 335,817-824), Appendix 00950 PURE FUNCTION kdL_DM(i,T) 00951 REAL(KIND=qprec),INTENT(IN)::T 00952 INTEGER,INTENT(IN):: i 00953 REAL(KIND=qprec)::kdL_DM 00954 00955 kdL_DM = 0 00956 SELECT CASE(i) 00957 CASE(H) 00958 kdL_DM = 4.69d-14*T**0.746d0*exp(-55065d0/T) 00959 END SELECT 00960 END FUNCTION kdL_DM 00961 00962 ! collisional dissociation rate for H2 + partner collision 00963 ! i = index for a particular collision partner 00964 ! H = 1 00965 ! H2 = 2 00966 ! see: Lepp & Schull (1983 ApJ) 00967 ! NOTE: these parameters are only good in low density limit n <~ 1e3 00968 ! for larger densities, these are lower bounds on the dissociation rate 00969 ! These lower density limit rates have been superceeded by newer work. 00970 ! This will probably not be used -- I have coded them here just so I can 00971 ! test or make comparisions with previous simulations done using these rates. 00972 PURE FUNCTION kdL_LS(i,T) 00973 REAL(KIND=qprec),INTENT(IN)::T 00974 INTEGER,INTENT(IN):: i 00975 REAL(KIND=qprec)::kdL_LS 00976 00977 kdL_LS = 0 00978 SELECT CASE(i) 00979 CASE(H) 00980 IF(T > 7390) THEN 00981 kdL_LS = 6.11d-14*exp(-2.93d4/T) 00982 ELSE 00983 kdL_LS = 2.67d-15*exp(-(6750d0/T)**2) 00984 END IF 00985 CASE(H2) 00986 IF(T > 7291) THEN 00987 kdL_LS = 5.22d-14*exp(-3.22d4/T) 00988 ELSE 00989 kdL_LS = 3.17d-15*exp(-4060.d0/T-(7500d0/T)**2) 00990 END IF 00991 END SELECT 00992 END FUNCTION kdL_LS 00993 END FUNCTION H2_diss 00994 00995 ! Recombination of H2 due to condensation on dust grains 00996 ! see: Hollenbach & McKee 1979 (1979ApJS...41..555H), page 564 00997 PURE FUNCTION H2_dust_Recomb(nHneuclei,nH,T,Tdust) !cm^-3*s^-1 00998 REAL(KIND=qPrec),INTENT(IN)::nHneuclei,nH,T,Tdust 00999 REAL(KIND=qPrec) H2_dust_Recomb 01000 REAL(KIND=qPrec) T2, Tdust2, fa 01001 01002 ! use fa = 0.5 as suggested in paper -- 01003 fa = 0.5 01004 T2 = T/100. 01005 Tdust2 = Tdust/100. 01006 01007 ! this function is only good to 10^5 01008 IF(T < 1.d5) THEN 01009 H2_dust_Recomb = nH*nHneuclei*(3.d-17*sqrt(T2)*fa/(1+0.4*sqrt(T2+Tdust2)+0.2*T2+0.08*T2**2)) 01010 ELSE 01011 H2_dust_recomb = 0. 01012 END IF 01013 01014 END FUNCTION H2_dust_Recomb 01015 01016 !Tablular form of function 01017 PURE FUNCTION H2_dust_Recomb_table(nHneuclei,nH,T) !cm^-3*s^-1 01018 REAL(KIND=qprec),INTENT(IN)::nHneuclei,nH,T 01019 REAL(KIND=qprec) H2_dust_Recomb_table,Tin,inc 01020 INTEGER :: lower,upper ! the lower,upper interpolation point 01021 01022 H2_dust_Recomb_table = zero 01023 IF(nHneuclei <=small .OR. nH <= small) RETURN 01024 Tin = log10(T) 01025 inc = (Tin-min_TAB)/delta_TAB 01026 lower = MAX(MIN(INT(inc),n_TAB),1) 01027 inc = inc-REAL(lower,qprec) 01028 upper = MAX(MIN(lower + 1,n_TAB),1) 01029 H2_dust_Recomb_table = nH*nHneuclei*(H2_dust_Recomb_TAB(lower)+(H2_dust_Recomb_TAB(upper)-H2_dust_Recomb_TAB(lower))*inc) 01030 END FUNCTION H2_dust_Recomb_table 01031 01032 ! Inonization / recombination rate calculations: 01033 ! General reference see: Marten & Szczerba 1997 (1997A&A...325.1132M) 01034 01035 ! radiative H+ + e- -> H recombination rate 01036 ! see: Mazzotta et. al. (1998A&AS..133..403M) 01037 ! see: Verner & Ferland (1996ApJS..103..467V) 01038 PURE FUNCTION H_recomb(nHII,ne,T) !cm^-3*s^-1 01039 REAL(KIND=qprec),INTENT(IN)::nHII,ne,T 01040 REAL(KIND=qprec) H_recomb,Tin 01041 REAL(KIND=qprec),PARAMETER :: a=7.982d-11, 01042 b=0.7480d0, 01043 To=3.148d0, 01044 T1=7.036d5 01045 REAL(KIND=qprec) :: c 01046 Tin=MAX(200.d0,T) 01047 c = sqrt(T/To) 01048 H_recomb = nHII*ne*(a*1.d0/(c*(1+c)**(1-b)*(1+sqrt(Tin/T1))**(1+b))) 01049 END FUNCTION H_recomb 01050 01051 !Tablular form of function 01052 PURE FUNCTION H_recomb_table(nHII,ne,T) !cm^-3*s^-1 01053 REAL(KIND=qprec),INTENT(IN)::nHII,ne,T 01054 REAL(KIND=qprec) H_recomb_table,Tin,inc 01055 INTEGER :: lower,upper ! the lower,upper interpolation point 01056 01057 H_recomb_table=zero 01058 IF(ne <= zero .OR. nHII <= zero .OR. T<200.0) RETURN 01059 01060 Tin = log10(T) 01061 inc = (Tin-min_TAB)/delta_TAB 01062 lower = MAX(MIN(INT(inc),n_TAB),1) 01063 inc = inc-REAL(lower,qprec) 01064 upper = MAX(MIN(lower + 1,n_TAB),1) 01065 H_recomb_table = ne*nHII*(H_recomb_TAB(lower)+(H_recomb_TAB(upper)-H_recomb_TAB(lower))*inc) 01066 ! H_recomb_table = 10**(H_recomb_TAB(lower)+(H_recomb_TAB(upper)-H_recomb_TAB(lower))*inc) 01067 ! H_recomb_table = ne*nHII*(H_recomb_table) 01068 END FUNCTION H_recomb_table 01069 01070 ! collisional ionization 01071 ! see: Mazzotta et. al. (1998A&AS..133..403M) 01072 ! see:Arnaud & Rothenflug (1985A&AS...60..425A) 01073 ! uses exponential integral function eone from specfun.f90 01074 FUNCTION H_ioniz(nH,ne,T) 01075 REAL(KIND=qprec),INTENT(IN)::nH,ne,T 01076 REAL(KIND=qprec) H_ioniz 01077 REAL(KIND=qprec),PARAMETER:: I=13.6, A=22.8, B=-12.0, C=1.9, D=-22.6, const=6.69d-7 01078 REAL(KIND=qprec),PARAMETER:: k_ev = 8.617342d-5 !boltzmann in ev/K 01079 REAL(KIND=qprec) x,f1,f2,F 01080 01081 ! H_ioniz=1.d-80 01082 H_ioniz=0.0 01083 IF(T<1000.) RETURN 01084 01085 x = I/(k_ev*T) 01086 f1=exp(x)*eone(x) 01087 !f1=df1(x) 01088 f2=dF2(x) 01089 01090 F=A*(1.-x*f1)+B*(1.+x-x*(2.+x)*f1)+C*f1+D*x*f2 01091 H_ioniz=nH*ne*(const*exp(-x)*F/((k_ev*T)**1.5d0*x)) 01092 END FUNCTION H_ioniz 01093 01094 !Tablular form of function 01095 PURE FUNCTION H_ioniz_table(nH,ne,T) !cm^-3*s^-1 01096 REAL(KIND=qprec),INTENT(IN)::nH,ne,T 01097 REAL(KIND=qprec) H_ioniz_table,Tin,inc 01098 INTEGER :: lower,upper ! the lower,upper interpolation point 01099 01100 H_ioniz_table=zero 01101 IF(ne <= zero .OR. nH <= zero .OR. T<1000.) RETURN 01102 01103 Tin = log10(T) 01104 inc = (Tin-min_TAB)/delta_TAB 01105 lower = MAX(MIN(INT(inc),n_TAB),1) 01106 inc = inc-REAL(lower,qprec) 01107 upper = MAX(MIN(lower + 1,n_TAB),1) 01108 H_ioniz_table = nH*ne*(H_ioniz_TAB(lower)+(H_ioniz_TAB(upper)-H_ioniz_TAB(lower))*inc) 01109 ! H_ioniz_table = 10**(H_ioniz_TAB(lower)+(H_ioniz_TAB(upper)-H_ioniz_TAB(lower))*inc) 01110 ! H_ioniz_table = nH*ne*(H_ioniz_table) 01111 END FUNCTION H_ioniz_table 01112 01113 ! collisional ionization 01114 ! see: Mazzotta et. al. (1998A&AS..133..403M) 01115 ! see:Arnaud & Rothenflug (1985A&AS...60..425A) 01116 ! uses exponential integral function eone from specfun.f90 01117 FUNCTION He_ioniz(nHe,ne,T) 01118 REAL(KIND=qprec),INTENT(IN)::nHe,ne,T 01119 REAL(KIND=qprec) He_ioniz 01120 REAL(KIND=qprec),PARAMETER::I=24.6, A=17.8,B=-11.0,C=7.0,D=-23.2, const=6.69d-7 01121 REAL(KIND=qprec),PARAMETER:: k_ev = 8.617342d-5 !boltzmann in ev/K 01122 REAL(KIND=qprec) x,f1,f2,F 01123 01124 IF(T>1000.) THEN 01125 x = I/(k_ev*T) 01126 f1=exp(x)*eone(x) 01127 !f1=df1(x) 01128 f2=dF2(x) 01129 01130 F=A*(1.-x*f1)+B*(1.+x-x*(2.+x)*f1)+C*f1+D*x*f2 01131 He_ioniz=nHe*ne*(const*exp(-x)*F/((k_ev*T)**1.5d0*x)) 01132 ELSE 01133 ! He_ioniz=1.d-134 01134 He_ioniz=0.0 01135 END IF 01136 END FUNCTION He_ioniz 01137 01138 !Tablular form of function 01139 PURE FUNCTION He_ioniz_table(nHeII,ne,T) !cm^-3*s^-1 01140 REAL(KIND=qprec),INTENT(IN)::nHeII,ne,T 01141 REAL(KIND=qprec) He_ioniz_table,Tin,inc 01142 INTEGER :: lower,upper ! the lower,upper interpolation point 01143 01144 01145 He_ioniz_table=zero 01146 IF(ne <= zero .OR. nHeII <= zero .OR. T<1000.) RETURN 01147 01148 Tin = log10(T) 01149 inc = (Tin-min_TAB)/delta_TAB 01150 lower = MAX(MIN(INT(inc),n_TAB),1) 01151 inc = inc-REAL(lower,qprec) 01152 upper = MAX(MIN(lower + 1,n_TAB),1) 01153 He_ioniz_table = ne*nHeII*(He_ioniz_TAB(lower)+(He_ioniz_TAB(upper)-He_ioniz_TAB(lower))*inc) 01154 ! He_ioniz_table = 10**(He_ioniz_TAB(lower)+(He_ioniz_TAB(upper)-He_ioniz_TAB(lower))*inc) 01155 ! He_ioniz_table = ne*nHeII*(He_ioniz_table) 01156 END FUNCTION He_ioniz_table 01157 01158 ! radiative He+ + e- -> He recombination rate 01159 ! do not forget ionization energy source term 01160 PURE FUNCTION He_recomb(nHeII, ne,T) !cm^-3*s^-1 01161 REAL(KIND=qprec),INTENT(IN):: nHeII, ne, T 01162 REAL(KIND=qprec) He_recomb,Tin 01163 Tin=MAX(200.d0,T) 01164 He_recomb = nHeII*ne*(He_rad_recomb(Tin) + He_dielec_recomb(Tin)) 01165 CONTAINS 01166 01167 ! radiative recombination 01168 ! see: Mazzotta et. al. (1998A&AS..133..403M) 01169 ! see: Verner & Ferland (1996ApJS..103..467V) 01170 PURE FUNCTION He_rad_recomb(T) 01171 REAL(KIND=qprec),INTENT(IN):: T 01172 REAL(KIND=qprec) He_rad_recomb 01173 REAL(KIND=qprec),DIMENSION(2),PARAMETER :: a=(/3.294d-11,9.356d-11/), 01174 b=(/0.6910d0,0.7892d0/), 01175 To=(/1.554d1,4.266d-2/), 01176 T1=(/3.676d7,4.677d6/) 01177 REAL(KIND=qprec) :: c 01178 01179 IF(T < 1.d6) THEN 01180 c = sqrt(T/To(1)) 01181 He_rad_recomb = a(1)*1.d0/(c*(1+c)**(1-b(1))*(1+sqrt(T/T1(1)))**(1+b(1))) 01182 ELSE 01183 c = sqrt(T/To(2)) 01184 He_rad_recomb = a(2)*1.d0/(c*(1+c)**(1-b(2))*(1+sqrt(T/T1(2)))**(1+b(2))) 01185 END IF 01186 END FUNCTION He_rad_recomb 01187 01188 ! dielectric recombination rate 01189 ! see: Mazzotta et. al. 1998A&AS..133..403M 01190 PURE FUNCTION He_dielec_recomb(T) 01191 REAL(KIND=qprec),INTENT(IN) :: T 01192 REAL(KIND=qprec) He_dielec_recomb 01193 REAL(KIND=qprec),PARAMETER :: c=0.11200E-08 ,E=39.70 01194 01195 He_dielec_recomb = T**(-3./2.) * c * exp(-E/T) 01196 END FUNCTION He_dielec_recomb 01197 END FUNCTION He_recomb 01198 01199 !Tablular form of function 01200 PURE FUNCTION He_recomb_table(nHeII,ne,T) !cm^-3*s^-1 01201 REAL(KIND=qprec),INTENT(IN)::nHeII,ne,T 01202 REAL(KIND=qprec) He_recomb_table,Tin,inc 01203 INTEGER :: lower,upper ! the lower,upper interpolation point 01204 01205 Tin = log10(T) 01206 inc = (Tin-min_TAB)/delta_TAB 01207 lower = MAX(MIN(INT(inc),n_TAB),1) 01208 inc = inc-REAL(lower,qprec) 01209 upper = MAX(MIN(lower + 1,n_TAB),1) 01210 He_recomb_table = ne*nHeII*(He_recomb_TAB(lower)+(He_recomb_TAB(upper)-He_recomb_TAB(lower))*inc) 01211 ! He_recomb_table = 10**(He_recomb_TAB(lower)+(He_recomb_TAB(upper)-He_recomb_TAB(lower))*inc) 01212 ! He_recomb_table = ne*nHeII*(He_recomb_table) 01213 END FUNCTION He_recomb_table 01214 01215 FUNCTION HeII_ioniz(nHeII,ne,T) 01216 REAL(KIND=qprec),INTENT(IN)::nHeII,ne,T 01217 REAL(KIND=qprec) HeII_ioniz 01218 REAL(KIND=qprec),PARAMETER::I=54.4, A=14.4,B=-5.6,C=1.9,D=-13.3, const=6.69d-7 01219 REAL(KIND=qprec),PARAMETER:: k_ev = 8.617342d-5 !boltzmann in ev/K 01220 REAL(KIND=qprec) x,f1,f2,F 01221 01222 01223 IF(T>1000.) THEN 01224 x = I/(k_ev*T) 01225 f1=exp(x)*eone(x) 01226 !f1=df1(x) 01227 f2=dF2(x) 01228 01229 F=A*(1.-x*f1)+B*(1.+x-x*(2.+x)*f1)+C*f1+D*x*f2 01230 HeII_ioniz=nHeII*ne*(const*exp(-x)*F/((k_ev*T)**1.5d0*x)) 01231 ELSE 01232 HeII_ioniz=1.d-80 01233 END IF 01234 END FUNCTION HeII_ioniz 01235 01236 PURE FUNCTION HeII_ioniz_table(nHeII,ne,T) !cm^-3*s^-1 01237 REAL(KIND=qprec),INTENT(IN)::nHeII,ne,T 01238 REAL(KIND=qprec) HeII_ioniz_table,Tin,inc 01239 INTEGER :: lower,upper ! the lower,upper interpolation point 01240 01241 Tin = log10(T) 01242 inc = (Tin-min_TAB)/delta_TAB 01243 lower = MAX(MIN(INT(inc),n_TAB),1) 01244 inc = inc-REAL(lower,qprec) 01245 upper = MAX(MIN(lower + 1,n_TAB),1) 01246 ! HeII_ioniz_table = ne*nHeII*(HeII_ioniz_TAB(lower)+(HeII_ioniz_TAB(upper)-HeII_ioniz_TAB(lower))*inc) 01247 HeII_ioniz_table = 10**(HeII_ioniz_TAB(lower)+(HeII_ioniz_TAB(upper)-HeII_ioniz_TAB(lower))*inc) 01248 HeII_ioniz_table = ne*nHeII*(HeII_ioniz_table) 01249 END FUNCTION HeII_ioniz_table 01250 01251 PURE FUNCTION HeIII_recomb(nHeIII,ne,T) !cm^-3*s^-1 01252 REAL(KIND=qprec),INTENT(IN)::nHeIII,ne,T 01253 REAL(KIND=qprec) HeIII_recomb,Tin 01254 REAL(KIND=qprec),PARAMETER :: a=1.891d-10, 01255 b=0.7524, 01256 To=9.370, 01257 T1=2.774d6 01258 REAL(KIND=qprec) :: c 01259 Tin=MAX(200.d0,T) 01260 c = sqrt(T/To) 01261 HeIII_recomb = nHeIII*ne*(a*1.d0/(c*(1+c)**(1-b)*(1+sqrt(Tin/T1))**(1+b))) 01262 END FUNCTION HeIII_recomb 01263 01264 PURE FUNCTION HeIII_recomb_table(nHeIII,ne,T) !cm^-3*s^-1 01265 REAL(KIND=qprec),INTENT(IN)::nHeIII,ne,T 01266 REAL(KIND=qprec) HeIII_recomb_table,Tin,inc 01267 INTEGER :: lower,upper ! the lower,upper interpolation point 01268 01269 Tin = log10(T) 01270 inc = (Tin-min_TAB)/delta_TAB 01271 lower = MAX(MIN(INT(inc),n_TAB),1) 01272 inc = inc-REAL(lower,qprec) 01273 upper = MAX(MIN(lower + 1,n_TAB),1) 01274 ! HeIII_recomb_table =ne*nHeIII*(HeIII_recomb_TAB(lower)+(HeIII_recomb_TAB(upper)-HeIII_recomb_TAB(lower))*inc) 01275 HeIII_recomb_table =10**(HeIII_recomb_TAB(lower)+(HeIII_recomb_TAB(upper)-HeIII_recomb_TAB(lower))*inc) 01276 HeIII_recomb_table =ne*nHeIII*(HeIII_recomb_table) 01277 END FUNCTION HeIII_recomb_table 01278 01279 SUBROUTINE zbar_calc(T,ion_fract,zbar) 01280 REAL(kind=qPrec), PARAMETER :: gHI=2.0,gHII=1.0,gOI=4.0,gOII=9.0,gNI=9.0,gNII=4.0 01281 REAL(kind=qPrec), PARAMETER :: En=14.53,Eh=13.6,Eo=13.62,ev_to_erg=1.602e-12,k_ev=8.617343e-5 01282 01283 REAL(kind=qPrec) :: OItoOII,NItoNII,g1O,g2O,g1N,g2N,c_O,c_N 01284 REAL(kind=qPrec),Intent(IN) :: T,ion_fract 01285 REAL(kind=qPrec) :: zbar 01286 REAL(kind=qPrec),Allocatable,Dimension(:) :: z_part,cz 01287 INTEGER :: num_metals,I 01288 01289 01290 num_metals = size(metal_abun) 01291 ALLOCATE(z_part(size(metal_abun))) 01292 ALLOCATE(cz(size(metal_abun))) 01293 cz = metal_states 01294 01295 !!!! Charge exchange 01296 g1O = gHII*gOI 01297 g2O = gHI*gOII 01298 g1N = gHII*gNI 01299 g2N = gHI*gNII 01300 OItoOII = (g1O/g2O)*(exp(-(Eh-Eo)/(k_ev*T)))*(((1.0-ion_fract)/ion_fract)) 01301 NItoNII = (g1N/g2N)*(exp(-(En-Eo)/(k_ev*T)))*(((1.0-ion_fract)/ion_fract)) 01302 c_O = 1.d0 / (1+OItoOII) 01303 c_N = 1.d0 / (1+NItoNII) 01304 !!!!!! 01305 01306 cz(2) = c_N ! Nitrogen will have both NI and NII in a ratio set by charge exhange 01307 cz(3) = c_O ! Oxygen likewise with OI and OIi 01308 01309 01310 01311 DO I = 1,num_metals 01312 z_part(I) = metal_abun(I)*metal_states(I)*cz(I) 01313 END DO 01314 zbar = SUM(z_part(1:num_metals)) / SUM(metal_abun(1:num_metals)) 01315 DEALLOCATE(z_part) 01316 DEALLOCATE(cz) 01317 END SUBROUTINE zbar_calc 01318 01319 01320 !analytic approximation to f2 integral for inonization rates in Arnaud et. al. 01321 PURE FUNCTION df2(x) 01322 REAL(KIND=qprec),INTENT(IN)::x 01323 REAL(KIND=qprec) df2 01324 REAL(KIND=qprec),DIMENSION(0:13),PARAMETER::p=(/1.d0,2.1658d2,2.0336d4,1.0911d6,3.7114d7,8.3963d8,1.2889d10, 01325 1.3449d11,9.4002d11,4.2571d12,1.1743d13,1.7549d13,1.0806d13,4.9776d11/) 01326 REAL(KIND=qprec),DIMENSION(0:14),PARAMETER::q=(/1.d0,2.1958d2,2.0984d4,1.1517d6,4.0349d7,9.49d8,1.5345d10, 01327 1.7182d11,1.3249d12,6.9071d12,2.3531d13,4.9432d13,5.7760d13,3.0225d13,3.3641d12/) 01328 REAL(KIND=qprec) ptot,qtot 01329 REAL(KIND=qprec),DIMENSION(0:14)::xin 01330 INTEGER j 01331 01332 FORALL(j=0:14) 01333 xin(j)=x**(-j) 01334 END FORALL 01335 01336 ptot=DOT_PRODUCT(xin(0:13),p(0:13)) 01337 qtot=DOT_PRODUCT(xin(0:14),q(0:14)) 01338 01339 df2=ptot/(x**2*qtot) 01340 END FUNCTION df2 01341 01342 ! approximation to exponential integral funciton, should be faster than 01343 ! exact calculation from netlib.f 01344 ! doesnt seem to work too well, not used 01345 PURE FUNCTION df1(x) 01346 REAL(KIND=qprec),INTENT(IN)::x 01347 REAL(KIND=qprec) df1 01348 INTEGER,PARAMETER::nterms=15 01349 REAL(KIND=qprec),DIMENSION(nterms)::term 01350 INTEGER j 01351 01352 IF(x < 0.02) THEN 01353 df1=exp(x)*(-log(x)-0.5772+x) 01354 ELSE IF(x < 1.5) THEN 01355 df1=log((x+1.)/x)-(0.36+0.03*(x+0.01)**(-0.5))/(x-1)**2 01356 ELSE IF(x < 10.) THEN 01357 df1=log((x+1.)/x)-(0.36+0.03*(x+0.01)**(0.5))/(x-1)**2 01358 ELSE 01359 term(1)=1./x 01360 DO j=2,nterms 01361 term(j)=-j*term(j-1)/x 01362 END DO 01363 !df1=1/x*(1. - 1./x + 2./(x**2) - 6./(x**3) + 24./(x**4) - 5.*24./(x**5)) 01364 df1=SUM(term) 01365 END IF 01366 END FUNCTION df1 01367 01368 SUBROUTINE Cool_Derivatives(q,dqdt,Temp,nvec,iCooling) 01369 REAL(KIND=qPrec), DIMENSION(1:NrHydroVars),INTENT(IN) :: q 01370 REAL(KIND=qPrec), DIMENSION(1:NrHydroVars), INTENT(OUT):: dqdt 01371 INTEGER, OPTIONAL :: iCooling 01372 REAL(KIND=qPrec), OPTIONAL, INTENT(IN) :: Temp 01373 REAL(KIND=qPREC), DIMENSION(0:nSpeciesHi), OPTIONAL, INTENT(IN) :: nvec 01374 REAL(KIND=qPrec) :: hion,nneuclei,nHneuc,H2diss,hrec,heion,herec,dm,OIcool,heIIIrec,heIIion, 01375 dustcool,H2cool,H2heat,dustheat,H2term,qHeat,cooling,eth,coolfrac,gamma,CDH2,CD,BBC,hex,EqIon,ionfrac,Hex_loss,ion_recomb_loss 01376 REAL(KIND=qPrec) :: minT=10**min_TAB, maxT=10**max_TAB 01377 REAL(KIND=qPrec) :: nneuc,ne,nH,mu 01378 INTEGER, PARAMETER :: ZCool=5 01379 REAL(KIND=qPrec) :: logne, logX 01380 01381 dqdt = 0d0 01382 01383 IF(.NOT. lchemion) RETURN 01384 01385 nneuc = nvec(iH2) + nvec(iH) + nvec(iHII) + nvec(iHe) + nvec(iHeII) + nvec(iHeIII) 01386 ne = nvec(iHII) + nvec(iHeII) + 2d0*nvec(iHeIII) + nmine*(nneuc) 01387 nH = nvec(iH) + nvec(iHII) 01388 logne = LOG10(ne) 01389 01390 H2diss=0d0 01391 hion = 0d0 01392 hrec = 0d0 01393 heion = 0d0 01394 herec = 0d0 01395 dm = 0d0 01396 OIcool=0d0 01397 dustcool=0d0 01398 H2cool=0d0 01399 H2term=0d0 01400 hex=0d0 01401 BBC=0d0 01402 01403 01404 IF(lHII) THEN 01405 ! hion = H_ioniz(nvec(iH),ne,Temp) 01406 ! hrec = H_recomb(nvec(iHII),ne,Temp) 01407 hion = H_ioniz_table(nvec(iH),ne,Temp) 01408 hrec = H_recomb_table(nvec(iHII),ne,Temp) 01409 END IF 01410 01411 IF(nvec(iHe)+nvec(iHeII)+nvec(iHeIII) > small .AND. lHeII) THEN 01412 ! heion = He_ioniz(nvec(iHe),ne,Temp) 01413 ! herec = He_recomb(nvec(iHeII),ne,Temp) 01414 heion = He_ioniz_table(nvec(iHe),ne,Temp) 01415 herec = He_recomb_table(nvec(iHeII),ne,Temp) 01416 END IF 01417 01418 IF((nvec(iHeIII)+nvec(iHeII)) > small .AND. lHeIII) THEN 01419 !heion = He_ioniz(nvec(iHe),ne,Temp) 01420 !herec = He_recomb(nvec(iHeII),ne,Temp) 01421 heIIion = HeII_ioniz_table(nvec(iHeII),ne,Temp) 01422 heIIIrec = HeIII_recomb_table(nvec(iHeIII),ne,Temp) 01423 END IF 01424 01425 IF(lneqcooling) THEN 01426 IF (nH .GT. 0d0) THEN 01427 ionfrac = nvec(iHII) / nH !(nvec(iHII) + nvec(iH)) 01428 logX = LOG10(ionfrac) 01429 dm = DMCoolingRate(Temp) !DM Cooling only from metal 01430 dm = dm * ionfrac 01431 ELSE 01432 dm =0d0 01433 END IF 01434 hex = HexCoolingRate(Temp) !Cooling from H excitiation 01435 END IF 01436 01437 01438 Hex_loss = ne*nvec(iH)*hex ! Energy lost by hyrdogen excitation 01439 metal_loss = nH*nH*dm ! Energy lost by metal excitiation 01440 ion_recomb_loss = IonH*hion + (Boltzmann*Temp)*hrec + IonHe*heion +(Boltzmann*Temp)*herec !Energy lost by ionization and recombination 01441 01442 Tfloor = 1d3 01443 01444 01445 IF(Temp > Tfloor) THEN 01446 cooling = Hex_loss + metal_loss + ion_recomb_loss 01447 01448 ! If this is ZCooling, and values are within ranges of ZCoolingTab then ZCoolingRate replaces metal_loss 01449 IF(PRESENT(iCooling)) THEN 01450 IF(iCooling==ZCool) THEN 01451 IF(Temp<16.5d3 .AND. Temp>2d3 .AND. logne<7.8d0 .AND. logne>0d0 & 01452 .AND. logX<0d0 .AND. logX>-2d0) cooling = cooling - metal_loss 01453 END IF 01454 END IF 01455 ELSE 01456 cooling = 0d0 01457 END IF 01458 01459 ! Provide the inhomogenious part of the conservation laws to be solved 01460 IF(lH2) dqdt(iH2) = H2term*muH2 01461 IF(lHII) dqdt(iHII) = (hion-hrec) 01462 ! IF(lH) dqdt(iH) = (-(hion-hrec) - 2.d0*H2term)*muH 01463 IF(lH) dqdt(iH) = -(hion-hrec) 01464 IF(lHe) dqdt(iHe) = (herec-heion) 01465 IF(lHeII) dqdt(iHeII) = -(herec-heion) 01466 ! IF(lHeII) dqdt(iHeII) = (heion-herec+heIIIrec)*muHe 01467 IF(lHeIII) dqdt(iHeIII) = 0d0 !(heIIion-heIIIrec)*muHe 01468 01469 IF(lneqcooling) dqdt(iE) = cooling/pscale ! do energy scaling 01470 01471 ! do density scaling 01472 dqdt(nSpeciesLO:nSpeciesHI) = dqdt(nSpeciesLO:nSpeciesHI)/nScale 01473 01474 ! do time scaling 01475 dqdt = dqdt*timescale 01476 01477 END SUBROUTINE Cool_Derivatives 01478 END MODULE NEQCoolingSrc