10 use physicalproperties
20 subroutine equcond(keq,ystate,ngas,nfv,mor,eps,dcell,fstrut,temp)
21 real(dp),
intent(out) :: keq
22 real(dp),
dimension(:),
intent(in) :: ystate
23 integer,
intent(in) :: ngas
24 integer,
intent(in) :: nfv
25 integer,
intent(in) :: mor(:)
26 real(dp),
intent(in) :: temp
27 real(dp),
intent(in) :: eps
28 real(dp),
intent(in) :: dcell
29 real(dp),
intent(in) :: fstrut
31 real(dp),
dimension(ngas) :: kg,yg,cc
35 cpg(1)=oxyheatcapacity(temp)
36 cpg(2)=nitrheatcapacity(temp)
37 cpg(3)=cdheatcapacity(temp)
38 cpg(4)=cypheatcapacity(temp)
44 cc(k)=cc(k)+ystate(ngas*(i-1)+k)
52 kg(i)=gasconductivity(temp,i)
53 if (abs(yg(i))<1e-6_dp)
then 59 kgas = weightedaverage(kg,yg)
61 kgas = extwassiljewa(kg,yg,tc,pc,mg,temp,1._dp)
63 kgas = lindsaybromley(kg,yg,tb,cpg,mg,temp)
65 kgas = pandeyprajapati(kg,yg,tb,mg,temp)
67 stop
'unknown gas conductivity model' 81 end subroutine equcond
89 real(dp) function weightedaverage(k,yin)
result(kmix)
90 real(dp),
dimension(:),
intent(in) :: k
91 real(dp),
dimension(:),
intent(in) :: yin
93 real(dp),
dimension(:),
allocatable :: y
98 print*,
'Input molar fractions to weightedAverage cannot be negative.' 105 kmix=kmix+yin(i)*k(i)
107 end function weightedaverage
117 real(dp) function extwassiljewa(k,yin,Tc,pc,M,T,eps)
result(kmix)
118 real(dp),
dimension(:),
intent(in) :: k
119 real(dp),
dimension(:),
intent(in) :: yin
120 real(dp),
dimension(:),
intent(in) :: Tc
121 real(dp),
dimension(:),
intent(in) :: pc
122 real(dp),
dimension(:),
intent(in) :: M
123 real(dp),
intent(in) :: T
124 real(dp),
intent(in) :: eps
127 real(dp),
dimension(:),
allocatable :: gam,y
128 real(dp),
dimension(:,:),
allocatable :: ktr,A
130 allocate(y(n),gam(n),ktr(n,n),a(n,n))
132 if (minval(y)<0)
then 133 print*,
'Input molar fractions to extWassiljewa cannot be negative.' 138 gam=210*(tc*m**3/pc**4)**(1/6._dp)
141 ktr(i,j)=gam(j)*(exp(0.0464_dp*t/tc(i))-exp(-0.2412_dp*t/tc(i)))/&
142 gam(i)/(exp(0.0464_dp*t/tc(j))-exp(-0.2412_dp*t/tc(j)))
143 a(i,j)=eps*(1+sqrt(ktr(i,j))*(m(i)/m(j))**0.25_dp)**2/&
144 sqrt(8*(1+m(i)/m(j)))
153 kmix=kmix+y(i)*k(i)/x
155 end function extwassiljewa
164 real(dp) function lindsaybromley(k,yin,Tb,cp,M,T)
result(kmix)
165 real(dp),
dimension(:),
intent(in) :: &
171 real(dp),
intent(in) :: &
175 real(dp),
dimension(:),
allocatable :: &
180 real(dp),
dimension(:,:),
allocatable :: A
182 allocate(y(n),cv(n),s(n),gam(n),a(n,n))
184 if (minval(y)<0)
then 185 print*,
'Input molar fractions to lindsayBromley cannot be negative.' 197 x=k(i)/k(j)*cp(j)/cp(i)*(9-5/gam(j))/(9-5/gam(i))
199 0.25_dp*(1+sqrt(x*(m(j)/m(i))**0.75_dp*(t+s(i))/(t+s(j))))**2*&
200 (t+sqrt(s(i)*s(j)))/(t+s(j))
209 kmix=kmix+y(i)*k(i)/x
211 end function lindsaybromley
220 real(dp) function pandeyprajapati(k,yin,Tb,M,T)
result(kmix)
221 real(dp),
dimension(:),
intent(in) :: &
226 real(dp),
intent(in) :: &
230 real(dp),
dimension(:),
allocatable :: &
233 real(dp),
dimension(:,:),
allocatable :: A
235 allocate(y(n),s(n),a(n,n))
237 if (minval(y)<0)
then 238 print*,
'Input molar fractions to pandeyPrajapati cannot be negative.' 248 a(i,j)=0.25_dp*(1+sqrt(k(i)/k(j)*(m(i)/m(j))**0.25_dp*&
249 (t+s(i))/(t+s(j))))**2*(t+sqrt(s(i)*s(j)))/(t+s(j))
258 kmix=kmix+y(i)*k(i)/x
260 end function pandeyprajapati
262 end module conductivity