120 int main(
int argc,
char **argv )
122 FILE *finput,*foutput;
123 tdata *extratargets, *listdummy;
125 int have_switch=0, have_target=0, have_command=0, have_output=0;
127 char *switchfile, *targetfile, *command, *outputname;
128 int c, errflag=0, listflag=0, fd=-1;
129 char tempfilename[15] =
"";
130 static char list[5] =
"List";
133 extern int optind, optopt;
135 while((c=getopt(argc,argv,
":c:hlo:s:t:")) != -1) {
178 fprintf(stderr,
"Option -%c requires an operand\n", optopt);
182 fprintf(stderr,
"Unrecognised option: -%c\n", optopt);
187 fprintf(stderr,
"usage: crange [-c COMMAND] [-h] [-l] [-o FILE] [-s switch.ini] [-t target.ini] <task file>\n");
188 fprintf(stderr,
" -c COMMAND = Execute this one-line command instead of reading it from a file.\n");
189 fprintf(stderr,
" -h = Print this help message and exit.\n");
190 fprintf(stderr,
" -l = Print the built-in target table and exit.\n");
191 fprintf(stderr,
" -o FILE = Write to this file instead of standard output.\n");
192 fprintf(stderr,
" -s switch.ini = Override the default switch values by reading this file.\n");
193 fprintf(stderr,
" -t target.ini = Override the default target values by reading this file.\n");
194 fprintf(stderr,
" <task file> = A file containing a list of tasks for crange. Required unless a command is specified with -c.\n");
202 extratargets = (have_target) ?
init_target(targetfile) : NULL;
203 if(argc-optind >= 1) {
204 sscanf(argv[optind],
"%s",inputname);
205 finput=fopen(inputname,
"r");
207 fprintf(stderr,
"Error opening task file!\n");
210 }
else if (have_command) {
214 strcpy(tempfilename,
"/tmp/cr.XXXXXX");
215 if ((fd = mkstemp(tempfilename)) == -1 || (finput = fdopen(fd,
"w+")) == NULL) {
218 unlink(tempfilename);
220 fprintf(stderr,
"%s: %s\n", tempfilename, strerror(errno));
226 fprintf(finput,
"%s\n",command);
229 fprintf(stderr,
"No task file specified!\n");
233 foutput=fopen(outputname,
"w");
235 fprintf(stderr,
"Error opening output file!\n");
242 run_range( finput, foutput, sswitch, extratargets );
245 if (have_command) unlink(tempfilename);
246 if (have_target) free(extratargets);
269 gsl_complex Cm, previousterm, term, sumterm;
272 term=GSL_COMPLEX_ONE;
273 sumterm=GSL_COMPLEX_ONE;
277 Cm=gsl_complex_rect(dm-1.0,0.0);
278 term=gsl_complex_mul(previousterm,
280 gsl_complex_div(gsl_complex_add(a,Cm),
281 gsl_complex_add(b,Cm)),gsl_complex_div_real(z,dm)));
282 sumterm=gsl_complex_add(sumterm,term);
283 }
while( gsl_complex_abs(term) > 1.0e-6 && gsl_complex_abs(previousterm) > 1.0e-6 );
302 double x, y, r, fj, cterm;
303 double aterm1,aterm2,aterm3;
304 double lterm1,lterm2,lterm3;
307 static double coeff[6]={76.18009172947146,
311 0.1208650973866179e-2,
312 -0.5395239384953e-5};
321 r=sqrt((x+5.5)*(x+5.5)+y*y);
323 aterm2=(x+0.5)*atan2(y,(x+5.5))-y;
324 lterm1=(x+0.5)*log(r);
325 lterm2=-y*atan2(y,(x+5.5)) - (x+5.5) + 0.5*log(2.0*
M_PI);
327 denom=1.000000000190015;
330 cterm=coeff[j-1]/((x+fj)*(x+fj)+y*y);
335 aterm3=atan2(num,denom);
336 lterm3 = 0.5*log(num*num + denom*denom);
337 GSL_SET_COMPLEX(&result,lterm1+lterm2+lterm3,aterm1+aterm2+aterm3);
339 result=gsl_complex_sub(gsl_complex_rect(log(
M_PI),0.0),
340 gsl_complex_add(result,
341 gsl_complex_log(gsl_complex_sin(
342 gsl_complex_mul_real(z,
M_PI)))));
368 double z23, z1, g, b2, b;
374 z23 = exp((2.0/3.0)*log(z0));
378 ((2.045)+ 2.000*exp(-0.04369*z0))*
381 exp((0.2643)*log(e1))*
382 exp(-(0.4171)*log(z0))
385 }
else if (z2 == 6.0) {
387 ((2.584)+ 1.910*exp(-0.03958*z0))*
390 exp((0.2433)*log(e1))*
391 exp(-(0.3969)*log(z0))
396 ((1.164 + 0.2319*exp(-0.004302*z2))+ 1.658*exp(-0.05170*z0))*
398 -(8.144+0.9876*log(z2))*
399 exp((0.3140+0.01072*log(z2))*log(e1))*
400 exp(-(0.5218+0.02521*log(z2))*log(z0))
405 capA=1.16-z2*(1.91e-03 - 1.26e-05*z2);
406 capB=(1.18-z2*(7.5e-03 - 4.53e-05*z2))/
ALPHA;
414 z1=z0*(1.0-capA*exp(-capB*b/z23));
437 double djdx(
double e1,
double z0,
double I0,
double f0,
double K,
short sswitch,
tdata *target)
451 f1=0.3070722*z1*z1*z2/(b2*a2);
453 J=0.5*f1*(f2-b2-delt+K)*(f0/I0);
545 double dedx(
double e1,
double rel0,
double z0,
double a1,
short sswitch,
tdata *target )
547 static double fva[10]={0.33,0.078,0.03,0.014,0.0084,
548 0.0053,0.0035,0.0025,0.0019,0.0014};
555 double S,REL,f1,f2,f3,f4,f6,f7,f8,f9;
557 double Spa=0.0,dpa,ldpa,l0,Lpa0,Lpa0s,Lpa1,Lpa;
566 f1=0.3070722*z1*z1*z2/(b2*a1*a2);
570 cadj=1.0e-6*(target->
iadj)*(target->
iadj)*etam2*(0.422377
571 +etam2*(0.0304043-etam2*0.00038106))
572 +1.0e-9*(target->
iadj)*(target->
iadj)*(target->
iadj)*etam2*(3.858019
573 +etam2*(-0.1667989 + etam2*0.00157955));
589 0.5*(-log(1.0+2.0*((5.4858e-04)*g/a1)) - ((5.4858e-04)*g/a1)*b2/(g*g)) :
592 (
ALPHA/
M_PI)*b2*(6.0822 + log(2.0*g)*(log(2.0*g)*(2.4167 + 0.3333*log(2.0*g))-8.0314)) :
598 Lpa0=(19.0/9.0)*(log(g/4.0) - 11.0/6.0);
599 Lpa0s=(19.0/9.0)*log(183.0*exp(-1.0/3.0*log(z2))/(1.0 + 4.0*6.25470095193633*183.0*exp(-1.0/3.0*log(z2))/g));
601 +(28.0*l0/9.0 - 446.0/27.0)*2.0*ldpa/(M_PI*
M_PI) + 14.0*4.0*ldpa*ldpa/(9.0*M_PI*M_PI));
603 Spa=4.08803936906434e-06*(z1*z1/a1)*(z2*z2/a2)*(1.0 + 1.0/z2)*g*Lpa;
606 Bbr=log(1.0 + 2.0*g*0.179524783764566/(exp((1.0/3.0)*log(a1)) + exp((1.0/3.0)*log(a2)))/a1);
607 Sbr=5.21721169334564e-07*(z1*z1/a1)*(z1*z1/a1)*(z2*z2/a2)*g*Bbr;
612 v=b*g/(
ALPHA*sqrt(z2));
624 fv=fva[i]*exp(-2.0*log(v/10.0));
625 f4=1.0 + 2.0*z1*fv/(sqrt(z2));
630 S=f1*(f2*f4+f3+f6-(delt/2.0)+f8+f9) + Sbr + Spa;
636 REL=f1*(f2*f4+f3+f6-delt/2.0 - 0.5*f7 +f8);
659 double b,cbar,X,X0,X1;
663 cbar=2.0*log(target->
iadj/target->
pla)+1.0;
664 b=sqrt(1.0 - 1.0/(g*g));
667 cbar-=2.303*log10(target->
etad);
668 X1-=0.5*log10(target->
etad);
669 X0-=0.5*log10(target->
etad);
672 return(target->
d0*exp(4.6052*(X-X0)));
673 }
else if (X >= X0 && X < X1) {
674 return(4.6052*X + exp(log(target->
a) + target->
m*log(X1-X)) - cbar);
676 return( 4.6052*X - cbar );
697 if( g < 1.8 )
return(0.0);
698 cbar=2.0*log(target->
iadj/target->
pla)+1.0;
699 b=sqrt(1.0 - 1.0/(g*g));
701 if( target->
etad > 0 ){
702 y+=log(target->
etad);
705 y0 = ( cbar >= 13.804 ) ? 1.502*cbar-11.52 : 9.212;
708 if(cbar<12.25)y0=9.212;
709 if(cbar<11.5)y0=8.751;
710 if(cbar<11.0)y0=8.291;
711 if(cbar<10.5)y0=7.830;
712 if(cbar<10.0)y0=7.370;
715 if(target->
iadj>=100.0){
717 y0 = ( cbar >= 5.215 ) ? 1.502*cbar-6.909 : 0.9212;
720 y0 = ( cbar >= 3.681 ) ? 1.502*cbar-4.606 : 0.9212;
723 if( y < y0 )
return( 0.0 );
725 if( y > y1 )
return(y-cbar);
727 dy3=(y1-y0)*(y1-y0)*(y1-y0);
729 return(y-cbar+a*(y1-y)*(y1-y)*(y1-y));
768 double bma(
double z1,
double b )
771 double b2,y,y2,sumr,fn,fn2,f3,f5;
772 double lambda, theta0, cosx, st;
780 msum=(int)(10.0*y) + 1;
785 sumr+=(1.0/(fn2+y2) - 1.0/fn2)/fn;
792 f3=-y2*(1.202+sumr) +
relbloch(z1,b,lambda,theta0);
793 Cz1=gsl_complex_rect(0.5,-y);
794 Cz2=gsl_complex_rect(1.0,y);
812 f5=0.5*b2*y*((1.725+(0.52-2.0*st)*
M_PI*cosx)
815 +y*((-2.696+b2*(4.569-0.494*b2))
816 +y*(-1.170+b2*(0.222+1.254*b2))
838 double relbloch(
double z12,
double b1,
double lambda,
double theta0 )
840 gsl_complex Cl1, Cl2;
841 gsl_complex Cf1,Cf2,Cf3,Cf4,Cf5,Cf6,Cf7;
842 gsl_complex Ci,Cnu,Cth,Cabgl;
843 double ge = 0.5772157;
846 double sigma,es,bloch2;
849 g=1.0/sqrt(1.0-b1*b1);
850 abgl=
ALPHA/(b1*g*lambda);
852 Ci=gsl_complex_rect(0.0,1.0);
853 Cnu=gsl_complex_rect(1.0, 2.0*nu);
854 Cth=gsl_complex_rect(theta0/2.0,0.0);
855 Cabgl=gsl_complex_rect(abgl,0.0);
856 Cf1=gsl_complex_div(gsl_complex_rect(0.0,-1.0),Cnu);
857 Cf2=gsl_complex_pow(Cth,Cnu);
858 Cf3=gsl_complex_pow(Cabgl,Cnu);
859 Cf4=gsl_complex_exp(gsl_complex_rect(0.0,2.0*sigma));
860 Cf5=gsl_complex_div(gsl_complex_rect(1.0,0.0),Cnu);
861 Cf6=gsl_complex_sub(gsl_complex_rect(log(theta0/2.0),0.0),Cf5);
862 Cf7=gsl_complex_add(gsl_complex_rect(log(4.0/abgl)+ge-1.0,0.0),Cf5);
863 Cl1=gsl_complex_mul(Cf1,
864 gsl_complex_sub(gsl_complex_mul_real(Cf2,2.0),gsl_complex_mul(Cf3,Cf4)));
865 Cl2=gsl_complex_mul(Cf1,
867 gsl_complex_mul_real(
868 gsl_complex_mul(Cf2,Cf6),2.0),gsl_complex_mul(Cf3,
869 gsl_complex_mul(Cf4,Cf7))));
871 bloch2=(
M_PI/2.0)*b1*b1*nu*es*(4.0*nu*log(2.0)*GSL_REAL(Cl1)
872 + (nu*
M_PI-1.0)*GSL_IMAG(Cl1)
873 + 2.0*nu*GSL_REAL(Cl2));
897 double lindhard(
double zz,
double aa,
double bb,
short sswitch )
899 gsl_complex Cexir, Cexis, Cedr, Ceds, Cske, Cmske, Cgrgs;
900 gsl_complex Clamr,Clams;
901 gsl_complex Caar, Cbbr, Caas, Cbbs, Czzr;
902 gsl_complex Cpi,Cone;
904 double compton=3.05573356675e-3;
905 double a3,eta,gg,rho,prh;
906 double a0,b0,a1,gz,z1,an,bn,anm1,bnm1,anp1,nn,asum,bsum,signk;
907 double frgr, fsgs, grgs, figi, H;
908 double k=0.0,k0,l,sk;
910 double dkm1=0.0,dmk=0.0,sd2,sdm2,sdd;
911 double term1,term2,term3,sumterm;
915 Cpi=gsl_complex_rect(
M_PI,0.0);
916 Cone=GSL_COMPLEX_ONE;
919 gg=1.0/sqrt(1.0 - bb*bb);
928 if ((gg < 10.0/rho) || !(sswitch &
SSWITCH_NS)) {
932 max = n == 1 ? 3 : 2;
939 l= (k>0) ? k : -k-1.0;
940 Cske=gsl_complex_rect(sk+1.0,eta);
941 Cexir=gsl_complex_sqrt(gsl_complex_div(
942 gsl_complex_rect(k,-eta/gg),gsl_complex_rect(sk,-eta)));
943 Cedr=gsl_complex_mul(Cexir,
944 gsl_complex_exp(gsl_complex_rect(0.0,
946 if( sswitch & SSWITCH_NS ){
947 Cmske=gsl_complex_rect(-sk+1.0,eta);
948 Cexis=gsl_complex_sqrt(gsl_complex_div(gsl_complex_rect(k,-eta/gg),gsl_complex_rect(-sk,-eta)));
949 Ceds=gsl_complex_mul(Cexis,
950 gsl_complex_exp(gsl_complex_rect(0.0,
954 Cbbr=gsl_complex_rect(2.0*sk + 1.0,0.0);
955 Cbbs=gsl_complex_rect(-2.0*sk + 1.0,0.0);
956 Czzr=gsl_complex_rect(0.0,2.0*prh);
957 Clamr=gsl_complex_mul(gsl_complex_mul(Cexir,gsl_complex_exp(gsl_complex_rect(0.0,-prh))),
959 Clams=gsl_complex_mul(gsl_complex_mul(Cexis,gsl_complex_exp(gsl_complex_rect(0.0,-prh))),
961 grgs=GSL_IMAG(Clamr)/GSL_IMAG(Clams);
965 +2.0*sk*log(2.0*prh));
966 if(cos(GSL_IMAG(Cgrgs))<1.0)grgs*= -1.0;
967 if(fabs(grgs) > 1.0e-9) {
968 frgr=sqrt((gg-1.0)/(gg+1.0))*GSL_REAL(Clamr)/GSL_IMAG(Clamr);
969 fsgs=sqrt((gg-1.0)/(gg+1.0))*GSL_REAL(Clams)/GSL_IMAG(Clams);
970 gz=-1.0*signk*(rho*gg+1.5*
ALPHA*zz);
973 a0=(1.0 + 2.0*fabs(k))*b0/(rho-gz);
982 bn=((rho-gz)*an +
ALPHA*z1*anm1/2.0)/(2.0*nn+2.0*fabs(k)+1.0);
983 anp1=((gz+rho)*bn -
ALPHA*z1*bnm1/2.0)/(2.0*nn + 2.0);
990 }
while(fabs(anm1/asum) > 1e-6 && fabs(bnm1/bsum) > 1e-6 );
991 figi= (k>0) ? asum/bsum : bsum/asum;
992 H=((frgr-figi)/(figi-fsgs))*grgs;
998 Ceds=GSL_COMPLEX_ZERO;
1000 r=1.0 + H*H + 2.0*H*(GSL_REAL(Cedr)*GSL_REAL(Ceds) + GSL_IMAG(Cedr)*GSL_IMAG(Cedr));
1001 dk[i]=gsl_complex_arg(gsl_complex_add(Cedr,gsl_complex_mul_real(Ceds,H)));
1004 sdm2=sin(dk[2]-dk[1]);
1005 term1 = k0*(k0+1.0)*sdm2*sdm2/(eta*eta*(2.0*k0 + 1.0));
1007 sd2=sin(dk[0]-dkm1);
1008 term1+=k0*(k0-1.0)*sd2*sd2/(eta*eta*(2.0*k0 - 1.0));
1010 sdd=sin(dk[0]-dk[2]);
1011 term2 = k0*sdd*sdd/(eta*eta*(4.0*k0*k0 - 1.0));
1012 term3 = term1 - 1.0/k0;
1013 sumterm += term2 + term3;
1017 pct = (term2 + term3)/sumterm;
1020 sumterm = -log(prh) - 0.2;
1022 lls = sumterm + 0.5*bb*bb;
1046 }
else if (x < 1.0) {
1047 while ((fabs(t) > 0.0001) || (n < 10)) {
1049 s+=t/((double)(n*n));
1054 while ((fabs(t) > 0.0001) || (n < 10)) {
1056 s+=t/((double)(n*n));
1059 return(
M_PI*
M_PI/12.0 + 0.5*log(x)*log(x) + s);
1086 double range(
double e,
double z1,
double a1,
short sswitch,
tdata *target,
int *tno )
1090 double de2,dr,dedx1,dedx2,dedx3,dedx4,e0,e1,e2,e3,e4,dt,dtk;
1097 for (k=0;k<
MAXAB;k++) {
1098 if ((trange[k].z1 == z1) && (trange[k].a1 == a1) &&
1099 (trange[k].sswitch == sswitch) &&
1112 while (k < MAXAB && trange[k].target != NULL) k++;
1121 dtk = difftime(ct,trange[l].timestamp);
1130 trange[*tno].
z1 = z1;
1131 trange[*tno].
a1 = a1;
1132 trange[*tno].
sswitch = sswitch;
1133 trange[*tno].
target = target;
1144 e1=e0+1.33998104*de2;
1145 dedx1=
dedx(e1,rel,z1,a1,sswitch,target);
1146 e2=e0+1.86113631*de2;
1147 dedx2=
dedx(e2,rel,z1,a1,sswitch,target);
1148 e3=e0+0.13886369*de2;
1149 dedx3=
dedx(e3,rel,z1,a1,sswitch,target);
1150 e4=e0+0.66001869*de2;
1151 dedx4=
dedx(e4,rel,z1,a1,sswitch,target);
1152 dr=de2*(0.65214515/dedx1 + 0.34785485/dedx2 + 0.34785485/dedx3
1153 + 0.65214515/dedx4);
1154 trange[*tno].
range[i] = trange[*tno].
range[i-1] + dr;
1183 double qrange(
double e,
double z1,
double a1,
short sswitch,
tdata *target )
1186 double de2,dr,dedx1,dedx2,dedx3,dedx4,e1,e2,e3,e4;
1187 double ei=8.0,lef,lei,en=1.0,pen,rel = 0.0;
1192 ra[0]=
benton(ei,z1,a1,target);
1196 entries = (double) (
MAXE - 1);
1199 en=exp(
M_LN10*(lei + ((
double)i)*lef/entries));
1201 e1=pen+1.33998104*de2;
1202 dedx1=
dedx(e1,rel,z1,a1,sswitch,target);
1203 e2=pen+1.86113631*de2;
1204 dedx2=
dedx(e2,rel,z1,a1,sswitch,target);
1205 e3=pen+0.13886369*de2;
1206 dedx3=
dedx(e3,rel,z1,a1,sswitch,target);
1207 e4=pen+0.66001869*de2;
1208 dedx4=
dedx(e4,rel,z1,a1,sswitch,target);
1209 dr=de2*(0.65214515/dedx1 + 0.34785485/dedx2 + 0.34785485/dedx3
1210 + 0.65214515/dedx4);
1211 ra[i] = ra[i-1] + dr;
1215 }
else if (e > 1.0 && e <= ei) {
1216 return(
benton(e,z1,a1,target) );
1218 return( e*
benton(en,z1,a1,target)/en );
1251 double term,logt,logi,loglambda;
1254 double cz[4],join[4];
1255 static double amn[3][4][4] = {
1257 {-8.72500, 1.88000, 0.741900, 0.752000},
1258 { 0.83090, 0.11140, -0.528800, -0.555890},
1259 {-0.13396, -0.06481, 0.126423, 0.128431},
1260 { 0.01262, 0.00540, -0.009341, -0.009306}
1263 {-7.6604e-01, 2.5398e+00, -2.4598e-01, 0.0000e+00},
1264 { 7.3736e-02, -3.1200e-01, 1.1548e-01, 0.0000e+00},
1265 { 4.0556e-02, 1.8664e-02, -9.9661e-03, 0.0000e+00},
1266 { 0.0000e+00, 0.0000e+00, 0.0000e+00, 0.0000e+00}
1269 {-8.0155e+00, 1.8371e+00, 4.5233e-02, -5.9898e-03},
1270 { 3.6916e-01, -1.4520e-02, -9.5873e-04, -5.2315e-04},
1271 {-1.4307e-02, -3.0142e-02, 7.1303e-03, -3.3802e-04},
1272 { 3.4718e-03, 2.3603e-03, -6.8538e-04, 3.9405e-05}
1275 static double czmn[4][4] = {
1276 {-6.000e-05, 5.252e-02, 1.285e-01, 0.000e+00},
1277 {-1.850e-03, 7.355e-02, 7.171e-02, -2.723e-02},
1278 {-7.930e-02, 3.323e-01, -1.234e-01, 1.530e-02},
1279 { 2.200e-01, 0.000e+00, 0.000e+00, 0.000e+00}
1281 static double cjoin[2][7] = {
1282 { 0.94, 20.19, -84.08, 132.98, -30.77, -102.29, 64.03},
1283 {12.62, -51.96, 199.14, -367.09, 327.06, -108.57, 0.00}
1287 logi=log( target->
iadj );
1288 for (l = 0; l < 2; l++) {
1289 join[l] = cjoin[l][m=6];
1290 while (m > 0) join[l] = 0.001*target->
iadj*join[l] + cjoin[l][--m];
1292 for (l = 0; l < 3; l++) {
1294 for (m = 3; m >= 0; m-- ) {
1295 term=amn[l][m][n=3];
1296 while (n > 0) term=term*logt + amn[l][m][--n];
1301 loglambda+=log( (target->
a2)/(target->
z2) );
1302 prnglo[l]=exp(loglambda);
1303 if (l == 1) prnglo[l]*=1.0e-03;
1306 b=sqrt(1.0 -1.0/(g*g));
1308 for (m = 0; m < 4; m++) {
1310 for (n = 3; n >= 0; n--) {
1317 if ( e < join[0]) l = 0;
1318 if ( e > join[1]) l = 2;
1321 }
else if ( x > 0.2 && x <= 2.0 ) {
1323 }
else if ( x > 2.0 && x <= 3.0 ) {
1325 }
else if ( x > 3.0 ) {
1328 bzz=(31.8+3.86*exp((5.0/8.0)*logi))
1329 *( (target->
a2)/(target->
z2) )*1.0e-06*exp((8.0/3.0)*log( z1 ));
1330 return(( (a1/cr)/(z1*z1) )*(prnglo[l] + bzz*cz[n]));
1348 double renergy(
double e,
double r0,
double z1,
double a1,
short sswitch,
tdata *target )
1356 rr =
range(e,z1,a1,sswitch,target,&tno);
1358 if( r < 0.0 )
return(0.0);
1363 if(r > trange[tno].
range[0]){
1365 while( r > trange[tno].range[i] ) i++;
1367 (trange[tno].range[i]-trange[tno].range[i-1]));
1401 double red1,red2,z1,a1;
1408 icols=fscanf(finput,
"%s %lf %lf %lf %lf %s\n",
1409 task,&red1,&red2,&z1,&a1,abs);
1410 if(icols < 6)
break;
1413 if(icols==6 && strncmp( target->
name,
"Unknown",
NAMEWIDTH ) != 0){
1414 if(strncmp( task,
"r", 1 )==0){
1415 out=
range(red1,z1,a1,sswitch,target,&tno);
1416 }
else if(strncmp( task,
"e", 1 )==0){
1417 out=
renergy(red1,red2,z1,a1,sswitch,target);
1418 }
else if(strncmp( task,
"d", 1 )==0){
1419 out=
dedx(red1,red2,z1,a1,sswitch,target);
1420 }
else if(strncmp( task,
"j", 1 )==0){
1421 out=
djdx(red1,z1,2.0,0.05,3.04,sswitch,target);
1424 fprintf(foutput,
"%f\n",out);
1444 #ifdef HAVE_INIPARSER_H
1447 if (access(switchfile,R_OK) == 0) {
1449 d = iniparser_load( switchfile );
1450 if (iniparser_getboolean(d,
"switch:barkas",0)) sswitch |=
SSWITCH_BA;
1451 if (iniparser_getboolean(d,
"switch:shell",0)) sswitch |=
SSWITCH_SH;
1452 if (iniparser_getboolean(d,
"switch:leung",0)) sswitch |=
SSWITCH_LE;
1453 if (iniparser_getboolean(d,
"switch:new delta",1)) sswitch |=
SSWITCH_ND;
1454 if (iniparser_getboolean(d,
"switch:new electron capture",0)) sswitch |=
SSWITCH_EC;
1455 if (iniparser_getboolean(d,
"switch:finite nuclear size",1)) sswitch |=
SSWITCH_NS;
1456 if (iniparser_getboolean(d,
"switch:kinematic",0)) sswitch |=
SSWITCH_KI;
1457 if (iniparser_getboolean(d,
"switch:radiative",0)) sswitch |=
SSWITCH_RA;
1458 if (iniparser_getboolean(d,
"switch:pair",0)) sswitch |=
SSWITCH_PA;
1459 if (iniparser_getboolean(d,
"switch:bremsstrahlung",0)) sswitch |=
SSWITCH_BR;
1460 iniparser_freedict(d);
1462 fprintf(stderr,
"Could not read switch file: %s. Using default crange switches.\n",switchfile);
1485 #ifdef HAVE_INIPARSER_H
1491 if (access(targetfile,R_OK) == 0) {
1492 d = iniparser_load( targetfile );
1493 n = iniparser_getnsec(d);
1496 section=iniparser_getsecname(d,i);
1497 sprintf(key,
"%s:%s",section,
"name");
1499 strncpy(table[i].name,iniparser_getstring(d,key,
"Unknown"),
NAMEWIDTH);
1500 sprintf(key,
"%s:%s",section,
"z2");
1502 table[i].
z2 = iniparser_getdouble(d,key,0.0);
1503 sprintf(key,
"%s:%s",section,
"a2");
1504 table[i].
a2 = iniparser_getdouble(d,key,0.0);
1505 sprintf(key,
"%s:%s",section,
"iadj");
1506 table[i].
iadj = iniparser_getdouble(d,key,0.0);
1507 sprintf(key,
"%s:%s",section,
"rho");
1508 table[i].
rho = iniparser_getdouble(d,key,0.0);
1509 sprintf(key,
"%s:%s",section,
"pla");
1510 table[i].
pla = iniparser_getdouble(d,key,0.0);
1511 sprintf(key,
"%s:%s",section,
"etad");
1512 table[i].
etad = iniparser_getdouble(d,key,0.0);
1513 sprintf(key,
"%s:%s",section,
"bind");
1514 table[i].
bind = iniparser_getdouble(d,key,0.0);
1515 sprintf(key,
"%s:%s",section,
"x0");
1516 table[i].
X0 = iniparser_getdouble(d,key,0.0);
1517 sprintf(key,
"%s:%s",section,
"x1");
1518 table[i].
X1 = iniparser_getdouble(d,key,0.0);
1519 sprintf(key,
"%s:%s",section,
"a");
1520 table[i].
a = iniparser_getdouble(d,key,0.0);
1521 sprintf(key,
"%s:%s",section,
"m");
1522 table[i].
m = iniparser_getdouble(d,key,0.0);
1523 sprintf(key,
"%s:%s",section,
"d0");
1524 table[i].
d0 = iniparser_getdouble(d,key,0.0);
1526 iniparser_freedict(d);
1530 strncpy(table[n].name,
"Unknown",
NAMEWIDTH);
1532 fprintf(stderr,
"Could not read target file: %s. Using built-in crange targets.\n",targetfile);
1548 for (k=0;k<
MAXAB;k++) {
1549 trange[k].
z1 = trange[k].
a1 = 0.0;
1553 for (l=0;l<
MAXE;l++) trange[k].
range[l] = 0.0;
1574 double decades, entries;
1578 return(exp(
M_LN10*(LOGTENEMIN + ((
double)i)*decades/entries)));
1599 static tdata targets[] = {
1601 {
"H" , 1.000, 1.008, 19.2, 8.3748e-05, 2.6300e-01, 1.0, 1.3606e-02, 1.8639, 3.2718, 0.14092, 5.7273, 0.00 },
1602 {
"He" , 2.000, 4.003, 41.8, 1.6632e-04, 2.6300e-01, 1.0, 7.7872e-02, 2.2017, 3.6122, 0.13443, 5.8347, 0.00 },
1603 {
"C" , 6.000, 12.011, 78.0, 2.0000e+00, 2.8803e+01, 0.0, 1.0251e+00, -0.0351, 2.4860, 0.20240, 3.0036, 0.10 },
1604 {
"N" , 7.000, 14.007, 82.0, 1.1653e-03, 6.9500e-01, 1.0, 1.4782e+00, 1.7378, 4.1323, 0.15349, 3.2125, 0.00 },
1605 {
"O" , 8.000, 15.999, 95.0, 1.3315e-03, 7.4400e-01, 0.0, 2.0359e+00, 1.7541, 4.3213, 0.11778, 3.2913, 0.00 },
1606 {
"Na" , 11.000, 22.990, 149.0, 9.7100e-01, 1.9641e+01, 0.0, 4.4097e+00, 0.2880, 3.1962, 0.07772, 3.6452, 0.08 },
1607 {
"Al" , 13.000, 26.980, 166.0, 2.6989e+00, 3.2860e+01, 0.0, 6.5929e+00, 0.1708, 3.0127, 0.08024, 3.6345, 0.12 },
1608 {
"Si" , 14.000, 28.090, 173.0, 2.3300e+00, 3.1055e+01, 0.0, 7.8751e+00, 0.2014, 2.8715, 0.14921, 3.2546, 0.14 },
1609 {
"P" , 15.000, 30.974, 173.0, 2.2000e+00, 2.9743e+01, 0.0, 9.2905e+00, 0.1696, 2.7815, 0.23610, 2.9158, 0.14 },
1610 {
"Ar" , 18.000, 39.948, 188.0, 1.6620e-03, 7.8900e-01, 1.0, 1.4382e+01, 1.7635, 4.4855, 0.19714, 2.9618, 0.00 },
1611 {
"Fe" , 26.000, 55.847, 286.0, 7.8740e+00, 5.5172e+01, 0.0, 3.4582e+01, -0.0012, 3.1531, 0.14680, 2.9632, 0.12 },
1612 {
"Ni" , 28.000, 58.690, 311.0, 8.9020e+00, 5.9385e+01, 0.0, 4.1324e+01, -0.0566, 3.1851, 0.16496, 2.8430, 0.10 },
1613 {
"Cu" , 29.000, 63.550, 323.0, 8.9600e+00, 5.8270e+01, 0.0, 4.4973e+01, -0.0254, 3.2792, 0.14339, 2.9044, 0.08 },
1614 {
"Ge" , 32.000, 72.610, 350.0, 5.3230e+00, 4.4141e+01, 0.0, 5.7047e+01, 0.3376, 3.6096, 0.07188, 3.3306, 0.14 },
1615 {
"Ag" , 47.000, 107.870, 470.0, 1.0500e+01, 6.1635e+01, 0.0, 1.4451e+02, 0.0657, 3.1074, 0.24585, 2.6899, 0.14 },
1616 {
"Ba" , 56.000, 137.330, 491.0, 3.5000e+00, 3.4425e+01, 0.0, 2.2118e+02, 0.4190, 3.4547, 0.18268, 2.8906, 0.14 },
1617 {
"Os" , 76.000, 190.200, 746.0, 2.2570e+01, 8.6537e+01, 0.0, 4.6939e+02, 0.0891, 3.5414, 0.12751, 2.9608, 0.10 },
1618 {
"Pt" , 78.000, 195.090, 790.0, 2.1450e+01, 8.4389e+01, 0.0, 5.0101e+02, 0.1484, 3.6212, 0.11128, 3.0417, 0.12 },
1619 {
"Au" , 79.000, 196.970, 790.0, 1.9320e+01, 8.0215e+01, 0.0, 5.1732e+02, 0.2021, 3.6979, 0.09756, 3.1101, 0.14 },
1620 {
"Pb" , 82.000, 207.200, 823.0, 1.1350e+01, 6.1072e+01, 0.0, 5.6834e+02, 0.3776, 3.8073, 0.09359, 3.1608, 0.14 },
1621 {
"U" , 92.000, 238.030, 890.0, 1.5370e+01, 7.7986e+01, 0.0, 7.6220e+02, 0.2260, 3.3721, 0.19677, 2.8171, 0.14 },
1622 {
"Air" , 7.312, 14.667, 85.4, 1.2048e-03, 7.0700e-01, 1.0, 1.7150e+00, 1.7418, 4.2759, 0.10914, 3.3994, 0.00 },
1623 {
"ArCO2" , 15.867, 34.892, 174.7, 1.7000e-03, 8.0120e-01, 1.0, 1.7150e+00, 1.7418, 4.2759, 0.10914, 3.3994, 0.00 },
1624 {
"BC-408" , 3.381, 6.248, 62.8, 1.0320e+00, 2.1534e+01, 1.0, 4.9530e-01, 0.1769, 2.6747, 0.11442, 3.3762, 0.00 },
1625 {
"BP-1" , 14.795, 32.575, 242.5, 3.0000e+00, 3.3636e+01, 0.0, 2.7489e+01, 0.0843, 3.6297, 0.06445, 3.3655, 0.00 },
1626 {
"CH2" , 2.667, 4.676, 57.4, 9.4000e-01, 2.1099e+01, 0.0, 3.5078e-01, 0.1370, 2.5177, 0.12108, 3.4292, 0.00 },
1627 {
"CO2" , 7.333, 14.670, 85.0, 1.8421e-03, 8.7400e-01, 1.0, 1.6990e+00, 1.6294, 4.1825, 0.11768, 3.3227, 0.00 },
1628 {
"CR-39" , 3.946, 7.413, 73.1, 1.4000e+00, 2.4780e+01, 0.0, 7.2426e-01, 0.1562, 2.6507, 0.12679, 3.3076, 0.00 },
1629 {
"CsI" , 54.000, 129.905, 553.1, 4.5100e+00, 3.9455e+01, 0.0, 2.0258e+02, 0.0395, 3.3353, 0.25381, 2.6657, 0.00 },
1630 {
"Halo" , 1.077, 1.238, 19.2, 4.0880e-24, 5.4342e-11, 1.0, 1.8554e-02, 2.0000, 3.5000, 0.13500, 5.7500, 0.00 },
1631 {
"Hosta" , 4.250, 8.091, 72.3, 1.4000e+00, 2.4722e-01, 0.0, 7.7212e-01, 0.1606, 2.6255, 0.12860, 3.3288, 0.00 },
1632 {
"ISM" , 1.077, 1.238, 19.2, 2.0440e-24, 3.8426e-11, 1.0, 1.8554e-02, 2.0000, 3.5000, 0.13500, 5.7500, 0.00 },
1633 {
"Kapton" , 5.026, 9.803, 75.9, 1.4200e+00, 2.4586e+01, 0.0, 9.1859e-01, 0.1509, 2.5631, 0.15972, 3.1912, 0.00 },
1634 {
"Kevlar" , 4.000, 7.567, 71.7, 1.4500e+00, 2.5229e+01, 0.0, 9.1859e-01, 0.1509, 2.5631, 0.15972, 3.1912, 0.00 },
1635 {
"Lexan" , 4.061, 7.706, 73.1, 1.2040e+00, 2.2915e+01, 0.0, 6.8789e-01, 0.1606, 2.6255, 0.12860, 3.3288, 0.00 },
1636 {
"LH2" , 1.000, 1.008, 21.8, 6.0000e-02, 7.0310e+00, 0.0, 1.3606e-02, 0.4759, 1.9215, 0.13483, 5.6249, 0.00 },
1637 {
"Mesh" , 21.400, 44.200, 223.0, 2.1500e+00, 2.9400e+01, 0.0, 3.4582e+01, -0.0012, 3.1531, 0.14680, 2.9632, 0.12 },
1638 {
"Mylar" , 4.456, 8.735, 78.7, 1.4000e+00, 2.4595e+01, 0.0, 8.4108e-01, 0.1562, 2.6507, 0.12679, 3.3067, 0.00 },
1639 {
"SiO2" , 10.000, 20.029, 139.2, 2.3200e+00, 3.1014e+01, 0.0, 3.9822e+00, 0.1385, 3.0025, 0.08408, 3.5064, 0.00 },
1640 {
"Teflon" , 8.000, 16.669, 99.1, 2.2000e+00, 2.9609e+01, 0.0, 2.1465e+00, 0.1648, 2.7404, 0.10606, 3.4046, 0.00 },
1641 {
"Water" , 3.333, 6.005, 75.0, 1.0000e+00, 2.1469e+01, 0.0, 6.8770e-01, 0.2400, 2.8004, 0.09116, 3.4773, 0.00 },
1643 {
"Unknown" , 0.000, 0.000, 0.0, 0.0000e+00, 0.0000e+00, 0.0, 0.0000e+00, 0.0000, 0.0000, 0.00000, 0.0000, 0.00 }
1646 if (strncmp(
"List", target,
NAMEWIDTH) == 0) {
1647 while (strncmp(targets[k].name,
"Unknown",
NAMEWIDTH) != 0) {
1652 return(&targets[0]);
1654 if (extratargets != NULL) {
1655 while (strncmp(extratargets[k].name,
"Unknown",
NAMEWIDTH) != 0) {
1656 if (strncmp(extratargets[k].name, target,
NAMEWIDTH) == 0)
return(&extratargets[k]);
1661 while (strncmp(targets[k].name,
"Unknown",
NAMEWIDTH) != 0) {
1662 if (strncmp(targets[k].name, target,
NAMEWIDTH) == 0)
break;
1668 return(&targets[k]);
1679 printf(
"[%s]\n",target->
name);
1680 printf(
"name = %10s ; Target name\n", target->
name);
1681 printf(
"z2 = %10.3f ; Mean nuclear charge\n", target->
z2);
1682 printf(
"a2 = %10.3f ; Mean nuclear mass\n", target->
a2);
1683 printf(
"iadj = %10.1f ; Logarithmic mean ionization potential [eV]\n", target->
iadj);
1684 printf(
"rho = %10.4e ; Density [g cm^-3]\n", target->
rho);
1685 printf(
"pla = %10.4e ; Plasma frequency [eV]\n", target->
pla);
1686 printf(
"etad = %10.1f ; Ratio of density to density at STP for gasses (zero for everything else)\n", target->
etad);
1687 printf(
"bind = %10.4e ; Total electronic binding energy [eV]\n", target->
bind);
1688 printf(
"X0 = %10.4f ; Density effect turn-on value\n", target->
X0);
1689 printf(
"X1 = %10.4f ; Density effect asymptotic bound\n", target->
X1);
1690 printf(
"a = %10.5f ; Density effect interpolation parameter\n", target->
a);
1691 printf(
"m = %10.4f ; Density effect interpolation parameter\n", target->
m);
1692 printf(
"d0 = %10.2f ; Low energy density effect value (zero for everything but conductors)\n", target->
d0);
tdata * init_target(char *targetfile)
Read optional target data file.
short init_switch(char *switchfile)
Initializes the value of of the switch bit field.
void print_target(tdata *target)
Prints a target table entry in INI format.
tdata * find_target(char *target, tdata *extratargets)
Finds target data corresponding to a target name.
double delta(double g, tdata *target)
Computes the density effect.
double lindhard(double zz, double aa, double bb, short sswitch)
Compute the Lindhard-Sørensen correction.
double olddelta(double g, tdata *target)
Computes an obsolete version of the density effect.
void init_table(void)
Initialize range-energy tables.
int main(int argc, char **argv)
Main crange program.
double qrange(double e, double z1, double a1, short sswitch, tdata *target)
Computes total range by direct integration of dE/dx.
double energy_table(int i)
Returns the energy corresponding to a value in a range table.
gsl_complex complex_lngamma(gsl_complex z)
Complex logarithm of the Gamma function.
Structure containing target data.
double relbloch(double z12, double b1, double lambda, double theta0)
Compute the relativistic Bloch correction.
range_table trange[MAXAB]
The range-energy table.
double dedx(double e1, double rel0, double z0, double a1, short sswitch, tdata *target)
Computes dE/dx.
void run_range(FILE *finput, FILE *foutput, short sswitch, tdata *extratargets)
Parses and executes the task list.
double bma(double z1, double b)
Computes the Bloch, Mott and Ahlen corrections.
Structure to store range tables.
double renergy(double e, double r0, double z1, double a1, short sswitch, tdata *target)
Extract energies from range tables.
double effective_charge(double z0, double e1, double z2, short sswitch)
Computes effective projectile charge.
double benton(double e, double z1, double a1, tdata *target)
Computes ranges at low energies.
gsl_complex complex_hyperg(gsl_complex a, gsl_complex b, gsl_complex z)
Confluent hypergeometric function.
double range(double e, double z1, double a1, short sswitch, tdata *target, int *tno)
Computes total range given initial energy.
double Fbrems(double x)
Compute a mathematical function related to bremsstrahlung.
double djdx(double e1, double z0, double I0, double f0, double K, short sswitch, tdata *target)
Computes primary ionization.