Commit a3d653d1 authored by Pavel's avatar Pavel

Update code for ELOSS fortran subroutine

parent 840bb8e1
...@@ -4,320 +4,256 @@ ...@@ -4,320 +4,256 @@
!== File: ELOSS.f90 !== File: ELOSS.f90
!== Date: 08/12/2000 !== Date: 08/12/2000
!====================================================================== !======================================================================
! ELOSS.f90 ! ELOSS.f90
! !
! FUNCTIONS/SUBROUTINES exported from ELOSS.dll: ! FUNCTIONS/SUBROUTINES exported from ELOSS.dll:
! ELOSS - subroutine ! ELOSS - subroutine
!======================================================================
!== Date 2024-05-03
! Some reworking of code:
! !
subroutine ELOSS(NEL,ZEL,AEL,WEL,DEN,ZP,AP,NE,ETAB,RE,ZW,AW) ! * remove obsolete fortran instructions;
! * replace common block with module;
! * change precision to standard C double;
! * remove restriction for 5 elements in material.
module m_eloss
use iso_c_binding, only: c_int,c_double
real(c_double),allocatable,dimension(:),private:: adj,as,z,f,um
real(c_double),private:: oz,ro,azm,aj,zion,aion, zmed,amed
integer,private :: jc
contains
subroutine ELOSS(NEL,ZEL,AEL,WEL,DEN,ZP,AP,NE,ETAB,RE,ZW,AW)&
bind(C,name='eloss_')
! Expose subroutine ELOSS to users of this DLL
!
!DEC$ ATTRIBUTES DLLEXPORT::ELOSS
IMPLICIT REAL*8(A-H,O-Z) IMPLICIT none
INTERFACE
SUBROUTINE rde(E, R, R1, R2, IRE)
REAL*8 E, R, R1, R2
INTEGER IRE
END SUBROUTINE rde
END INTERFACE
! Variables ! Variables
! Input: ! Input:
INTEGER NEL, NE INTEGER(c_int),intent(in) ::NEL ! number of elements in the material
REAL*8 ZEL(NEL),AEL(NEL),WEL(NEL), DEN,ZP,AP, ETAB(NE) INTEGER(c_int) ::NE ! number of points in E -- Range table
real(c_double),dimension(nel),intent(in):: zel
real(c_double),dimension(nel),intent(in):: ael
real(c_double),dimension(nel),intent(in):: wel
real(c_double),intent(in) :: den
real(c_double),intent(in) :: ZP,AP
real(c_double),dimension(ne),intent(in):: Etab
! Output: ! Output:
REAL*8 RE(NE),ZW,AW real(c_double),dimension(ne),intent(out) :: RE
REAL(c_double),intent(out) :: ZW,AW
! Local
COMMON adj(5),as(5),z(5),f(5),um(5),oz,ro,azm,aj,zion,aion,&
zmed,amed,jc
!!$ local variables
integer :: i,j,ire
real(c_double) :: rel,rel1,uma,umed,zam,zme
! Body of ELOSS ! Body of ELOSS
zion = ZP zion = ZP
aion = AP aion = AP
jc = NEL jc = NEL
ro = DEN ro = DEN
allocate(adj(jc),as(jc),z(jc),f(jc),um(jc))
!== two parameters with question
ire = 1 !== two parameters with question
oz = 1. ire = 1
DO i = 1, NEL oz = 1.
as(i) = AEL(i)
z(i) = ZEL(i) as=ael;z = ZEL; um = WEL
um(i) = WEL(i) DO concurrent (i = 1:jc)
IF(z(i) .LE. 13.) THEN IF(z(i) .LE. 13.) THEN
adj(i) = 12.*z(i) + 7. adj(i) = 12.*z(i) + 7.
ELSE ELSE
adj(i) = 9.76*z(i) + 58.8 / z(i)**0.19 adj(i) = 9.76*z(i) + 58.8 / z(i)**0.19
END IF END IF
END DO END DO
uma=0. uma=sum(um*as)
zam=0. zme=sum(um*z)
umed=0. zam=zme/uma; azm=1./zam
zme=0. umed=sum(um)
aj=0. aj=exp(sum(um*z*log(adj))/uma*azm)
DO i = 1, jc
uma=uma+um(i)*as(i) amed=uma/umed
umed=umed+um(i) zmed=zme/umed
zme=zme+um(i)*z(i)
END DO f=um*as/uma
amed=uma/umed ZW = zmed
zmed=zme/umed AW = amed
DO i = 1, jc DO j = 1, NE
f(i)=um(i)*as(i)/uma !== result [mg/cm^2]
END DO CALL rde(ETAB(j), RE(j), rel1, rel, ire)
END DO
DO i = 1, jc deallocate(adj,as,z,f,um)
zam=zam+f(i)*z(i)/as(i)
END DO
azm=1./zam
aj=0.
DO i = 1, jc
aj=aj+f(i)*z(i)*log(adj(i))/as(i)
END DO
aj=aj*azm
aj=exp(aj)
ZW = zmed
AW = amed
DO j = 1, NE
!== result [mg/cm^2]
CALL rde(ETAB(j), RE(j), rel1, rel, ire)
END DO
end subroutine ELOSS end subroutine ELOSS
SUBROUTINE rde(e,range,rel1,rel,ix) SUBROUTINE rde(e,range,rel1,rel,ix)
IMPLICIT REAL*8(A-H,O-Z) ! calculates range and de/dx for compounds
IMPLICIT none !REAL(c_double) (A-H,O-Z)
INTERFACE
REAL*8 FUNCTION c(x) real(c_double),intent(in) ::E
REAL*8 x integer,intent(in) ::ix
END FUNCTION c real(c_double),intent(out) ::range
real(c_double),intent(out) ::rel1
REAL*8 FUNCTION c1(x) real(c_double),intent(out) ::rel
REAL*8 x
END FUNCTION c1 integer :: i,j,k,kk
END INTERFACE
real(c_double),parameter,dimension(3,3)::a=reshape([&
-0.75265, .073736, .040556,&
+2.53980,-.312000, .018664, &
! calculates range and de/dx for compounds -0.24598, .115480,-.0099661],shape=[3,3])
dimension ala(3),ala1(3) real(c_double),parameter,dimension(4,4)::a1=reshape([&
dimension a(3,3),a1(4,4),a2(4,4),b(2,3),cc(5) -8.0155 , 0.36916 ,-1.4307e-2, 3.4718e-3, &
common adj(5),as(5),z(5),f(5),um(5),oz,ro,azm,aj,zion,aion,& +1.8371 ,-1.4520e-2,-3.0142e-2, 2.3603e-3, &
zmed,amed,jc -4.5233e-2,-9.5873e-4, 7.1303e-3,-6.8538e-4, &
dimension alaj1(4),altau1(4) -5.9898e-3,-5.2315e-4,-3.3802e-4, 3.9405e-5],shape=[4,4])
real(c_double),parameter,dimension(4,4)::a2=reshape([&
data a/ -.75265, .073736, .040556, 2.5398,-.312,.018664, & -8.725000, 0.8309000,-0.13396000, 0.012625, &
-.24598, .11548, -.0099661/ -1.879700, 0.1113900,-0.06480800, 0.0054043, &
-0.741920,-0.5288050, 0.12642320,-9.341420e-3,&
data a1/-8.0155, 0.36916, -1.4307e-02, 3.4718e-03, & -0.752005,-0.5558937, 0.12843119,-9.306336e-3 ],shape=[4,4])
1.8371, -1.452e-02, -3.0142e-02, 2.3603e-03, & real(c_double),parameter,dimension(2,3)::b=reshape([&
0.045233, -9.5873e-04, 7.1303e-03, -6.8538e-04, & +4.223770e-07, 3.858019e-09, 3.04043e-08,&
-5.9898e-03,-5.2315e-04, -3.3802e-04, 3.9405e-05/ -1.667989e-10,-3.810600e-10, 1.57955e-12],shape=[2,3])
data a2/-8.725, 0.8309, -0.13396, 0.012625, & real(c_double),allocatable,dimension(:) :: cc
1.8797, 0.11139, -0.064808, 0.0054043, & real(c_double),dimension(3) :: ala, ala1
0.74192, -0.528805, 0.1264232, -0.00934142, & real(c_double) :: alaa,alaa1
0.752005, -0.5558937, 0.12843119, -0.009306336/ real(c_double),dimension(4) :: alaj1, altau1
real(c_double) :: abet, alaj, altau, ank, bep, beta, bz, bz1,cbet, t, tau, tt
data b/ 0.422377e-06, 3.858019e-09, 0.0304043e-06, -0.1667989e-09, & real(c_double),parameter :: coefa=3., coefb=1.
-0.00038106e-06, 0.00157955e-09/ real(c_double) :: alim1, alim2, del1, del2, en, hi, om, rel2, rel3, zef, zef2
if(e.le.0.002)then if(e.le.0.002)then
range=0. range=0.
rel1=0. rel1=0.
rel=0. rel=0.
return return
endif endif
! this a flag en = e/aion
ibis=-1 tau = en*1.008
en = e/aion altau = log(tau)
tau = en*1.008 alaj = log(aj)
DO concurrent (kk = 1:4)
altau = log(tau) alaj1(kk)=alaj**(kk-1)
alaj = log(aj) altau1(kk)=altau**(kk-1)
alaj1(1)= 1. END DO
altau1(1)= 1. t = tau/938.256
DO kk = 2,4 tt = 2.*t + t*t
alaj1(kk)=alaj**(kk-1) beta= sqrt(tt)/(1.+t)
altau1(kk)=altau**(kk-1)
END DO
t = tau/938.256 ala(1) = azm*exp(dot_product(alaj1,matmul(a2,altau1)))
tt = 2.*t + t*t ala(2) = azm*exp(dot_product(alaj1(1:3),matmul(a,altau1(1:3))))/1000.
beta= sqrt(tt)/(1.+t) ala(3) = azm*exp(dot_product(alaj1,matmul(a1,altau1)))
ala1(1)=ala(1)*dot_product(alaj1,matmul(a2(:,2:4),altau1(1:3)*[1,2,3]))/tau
s1=0. ala1(2)=ala(2)*dot_product(alaj1(1:3),matmul(a(:,2:3),altau1(1:2)*[1,2]))/tau
DO i=1,4 ala1(3)=ala(3)*dot_product(alaj1,matmul(a1(:,2:4),altau1(1:3)*[1,2,3]))/tau
DO j=1,4
s1 = s1 + a2(j,i)*alaj1(j)*altau1(i) associate(ca=>coefa*(.98 - en),cb=>coefb*(8.0 - en))
END DO associate(tca=>(1+tanh(ca))/2.,tcb=>(1+tanh(cb))/2.)
END DO alaa=(ala(1)*tca + ala(2)*(1.0-tca))*tcb + ala(3)*(1.0-tcb)
ala(1) = azm*exp(s1) alim1=0.; alim2=0.
if(-ca.lt.85)then
s2=0. alim1=1.008*(cosh(ca)**(-2))
DO i=2,4 endif
DO j=1,4 if(-cb.lt.85)then
s2 = s2 + (i-1)*a2(j,i)*alaj1(j)*altau1(i-1) alim2=1.008*(cosh(cb)**(-2))
END DO endif
END DO
ala1(1)=ala(1)*s2/tau alaa1=(ala1(1)*tca+ ala1(2)*(1.0-tca))*tcb + ala1(3)*(1.0-tcb)+&
coefa*tcb* (ala(2)-ala(1))*alim1/2.+ &
s1=0. coefb*(ala(3)-(ala(1)*tca+ala(2)*(1.-tca)))*alim2/2.
DO i=1,4 end associate
DO j=1,4 end associate
s1 = s1 + a1(j,i)*alaj1(j)*altau1(i)
END DO hi=137.*beta/zion
END DO bz=(31.8+3.86*(aj**.625))*azm* 1e-6 *(zion**2.666667)*c(hi)
bz1=(4.357+.5288*(aj**.625))*azm* .001 *(zion**1.666667)*c1(hi)
ala(3) = azm*exp(s1) bep=beta*beta
rel1=zion*zion/(alaa1+bz1*((1.- beta**2)**1.5)/931.141/beta)/1000.
s2=0. range=(alaa+bz)*aion/(1.008*zion**2)*1000.
DO i=2,4
DO j=1,4 ! Atention!! this version do not work correctly for ix.ne.1
s2 = s2 + (i-1)*a1(j,i)*alaj1(j)*altau1(i-1) if(ix.ne.1) then
END DO return
END DO end if
ala1(3)=ala(3)*s2/tau
allocate(cc(size(as)))
s1=0.
DO i=1,3 ank=.153574*ro/azm
DO j=1,3 ! z23=zion**.666667
s1 = s1 + a(j,i)*alaj1(j)*altau1(i) abet=beta*125.*(zion** (-2.0/3.0))
END DO zef=zion*(1.-exp(-abet))
END DO zef2=zef*zef
ala(2)=azm*exp(s1)/1000. om= 1022000. *bep/(1.-bep)
cbet=0.
s2=0. DO k=1,jc
DO i=2,3 cc(k)=0.
DO j=1,3 DO i=1,2
s2 = s2 + (i-1)*a(j,i)*alaj1(j)*altau1(i-1) DO j=1,3
END DO cc(k)=cc(k)+b(i,j)*((1./bep-1.)**j)*((adj(k)**(i+1)))
END DO END DO
ala1(2)=ala(2)*s2/tau END DO
cbet=cbet+f(k)*cc(k)/as(k)
25 continue END DO
coefa=3. cbet=cbet*azm
coefb=1. del1=ank*zef2*(log(om/oz)-bep)/bep/ro
del2=2.*ank*zef2*(log(oz/aj)-cbet)/bep/ro
alaa=(ala(1)*(tanh(coefa*(.98 - en))+1.)/2. & rel2=rel1-del1
+ ala(2)*(tanh(coefa*(en - .98))+1.)/2.) & rel3=rel1-rel2+del2
* (tanh(coefb*(8.0 - en))+1.)/2. & !!$ here is a replasement for an Arithmetic IF construct
+ ala(3)*(tanh(coefb*(en - 8.))+1.)/2. if (del1 <=0.0 ) then
rel=rel1
else if (del1+del2-rel2<=0) then
alim1=0. rel=rel3
alim2=0. else if(del1.lt.rel1)then
if(coefa*(en-.98).lt.85)then rel=rel2
alim1=1.008/cosh(coefa*(.98-en))/cosh(coefa*(.98-en)) else
endif rel=rel1
if(coefb*(en-8.).lt.85)then endif
alim2=1.008/cosh(coefb*(8.-en))/cosh(coefb*(8.-en)) deallocate(cc)
endif return
END SUBROUTINE rde
alaa1=(ala1(1)*(tanh(coefa*(.98-en))+1.)/2.+ &
ala1(2)*(tanh(coefa*(en-.98))+1.)/2.)* &
(tanh(coefb*(8.-en))+1.)/2.+ & pure FUNCTION c(x) result(res)
ala1(3)*(tanh(coefb*(en-8.))+1.)/2.+ & REAL(c_double),intent(in) :: x
coefa/2.*(tanh(coefb*(8.-en))+1.)/2.* & real(c_double) :: res
(ala(2)*alim1-ala(1)*alim1)+ & IF(x .LE. 0.2) THEN
coefb/2.*(ala(3)*alim2- & res = -0.00006 + (0.05252 + 0.12857*x)*x
(ala(1)*(tanh(coefa*(.98-en))+1.)/2.+ & ELSE IF(x .LE. 2.) THEN
ala(2)*(tanh(coefa*(en-.98))+1.)/2.)*alim2) res = -0.00185 + (0.07355 + (0.07172 - 0.02723*x)*x)*x
ELSE IF(x .LE. 3.) THEN
res = -0.0793 + (0.3323 - (0.1234 - 0.0153*x)*x)*x
hi=137.*beta/zion ELSE
bz=(31.8+3.86*(aj**.625))*azm*.000001 res = 0.22
bz=bz*(zion**2.666667)*c(hi) END IF
bz1=(4.357+.5288*(aj**.625))*azm*.001 END FUNCTION c
bz1=bz1*(zion**1.666667)*c1(hi)
bep=beta*beta pure FUNCTION c1(x) result(res)
rel1=zion*zion/(alaa1+bz1*((1.-bep)**1.5)/931.141/beta) REAL(c_double),intent(in) :: x
rel1=rel1/1000. real(c_double) :: res
range=(alaa+bz)*aion/1.008/zion/zion IF(x .LE. 0.2) THEN
range=range*1000. res = 0.05252+.25714*x
ELSE IF(x .LE. 2.0) THEN
! Atention!! this version do not work correctly for ix.ne.1 res = 0.07355 + (0.14344 - 0.08169*x)*x
if(ix.ne.1)return ELSE IF(x .LE. 3.0) THEN
ank=.153574*ro/azm res = 0.3323 - (0.2468 - 0.0459*x)*x
z23=zion**.666667 ELSE
abet=beta*125./z23 res = 0.
zef=zion*(1.-exp(-abet)) END IF
zef2=zef*zef
om=1022000.*bep/(1.-bep) END FUNCTION c1
cbet=0.
DO k=1,jc end module m_eloss
cc(k)=0.
DO i=1,2
DO j=1,3
cc(k)=cc(k)+b(i,j)*((1./bep-1.)**j)*((adj(k)**(i+1)))
END DO
END DO
cbet=cbet+f(k)*cc(k)/as(k)
END DO
cbet=cbet*azm
del1=ank*zef2*(log(om/oz)-bep)/bep
del1=del1/ro
del2=2.*ank*zef2*(log(oz/aj)-cbet)/bep
del2=del2/ro
rel2=rel1-del1
rel3=rel1-rel2+del2
if(del1)8,8,9
8 rel=rel1
goto 12
9 if(del1+del2-rel2)10,10,11
10 rel=rel3
goto 12
11 if(del1.lt.rel1)then
rel=rel2
else
rel=rel1
endif
12 return
END
REAL*8 FUNCTION c(x)
REAL*8 x
IF(x .LE. 0.2) THEN
c = -0.00006 + (0.05252 + 0.12857*x)*x
ELSE IF(x .LE. 2.) THEN
c = -0.00185 + (0.07355 + (0.07172 - 0.02723*x)*x)*x
ELSE IF(x .LE. 3.) THEN
c = -0.0793 + (0.3323 - (0.1234 - 0.0153*x)*x)*x
ELSE
c = 0.22
END IF
END
REAL*8 FUNCTION c1(x)
REAL*8 x
IF(x .LE. 0.2) THEN
c1 = 0.05252+.25714*x
ELSE IF(x .LE. 2.0) THEN
c1 = 0.07355 + (0.14344 - 0.08169*x)*x
ELSE IF(x .LE. 3.0) THEN
c1 = 0.3323 - (0.2468 - 0.0459*x)*x
ELSE
c1 = 0.
END IF
END
Markdown is supported
0% or
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment