17 implicit real *8 (a-h,o-z)
19 include
'common_phs.cmn'
23 character*255 infile,outdir
25 character*80 xname,outz,tmat,pf,dnr
27 character*2 ciband,crh,cisd,cset
32 open(21,file=
'phs.dir',status=
'old')
36 len1=
index(infile,
' ')-1
37 len2=
index(outdir,
' ')-1
39 open(5, file=infile(1:len1),status=
'old')
45 write(*,*) ilam1,ilam2
51 read(5,28)jfunc(i),kfunc(i)
59 read(5,30)(epar(i,ik),ik=1,6)
62 read(5,30)xww(i,j),xn1(i,j,1),xn2(i,j,1),xdx(i,j,1)
66 read(5,28)jrgm(i),irh(i),iset(i)
72 read(5,30)xww(i,j),xn1(i,j,k),xn2(i,j,k),xrg(i,j,k),
73 1 xsig(i,j,k),xnpar(i,j,k),xdx(i,j,k)
86 if(irh(isd) .ne. khum)
then
94 ixwwz=dnint(xww(isd,il)*1.0e3+0.01)
98 call convtc(irh(isd),2,crh)
99 call convtc(iset(isd),2,cset)
100 call convtc(ixwwz,4,cxwwz)
102 write(*,*) ixwwz, xww(isd,il), cxwwz,
' ',isd,
' of ',isd2
104 xname=
'wl'//cxwwz//
'sd'//cisd//
'rh'//crh//
'_set'//cset//
'.dat'
105 kx1=
index(xname,
' ')-1
107 outz=outdir(1:len2)//
'/phs_o'//xname(1:kx1)
108 tmat=outdir(1:len2)//
'/phs_t'//xname(1:kx1)
109 pf=outdir(1:len2)//
'/phs_p'//xname(1:kx1)
110 dnr=outdir(1:len2)//
'/phs_dn'//xname(1:kx1)
112 ko1=
index(outz,
' ')-1
113 kt1=
index(tmat,
' ')-1
118 open(unit=7,file=tmat(1:kt1),status=
'unknown')
119 open(unit=8,file=pf(1:kp1),status=
'unknown')
120 open(unit=4,file=dnr(1:kn1),status=
'unknown')
140 36
format(i5,1pe10.3)
141 40
format(10x,1p7e10.3)
142 50
format(
'jfunc is=',i2,
' max. value allowed is 5')
147 subroutine flbfr(il,isd,ifunc)
151 implicit real *8 (a-h,o-z)
153 include
'common_phs.cmn'
163 if(ifunc.ge.3 .and. ifunc.le.5)
then
167 ddd(mm+1)=xrg(isd,il,i)
168 ddd(mm+2)=xsig(isd,il,i)
169 ddd(mm+3)=xnpar(isd,il,i)
178 subroutine pfunc(il,isd)
182 implicit real *8 (a-h,o-z)
184 include
'common_phs.cmn'
200 call flbfr(il,isd,ifunc)
202 call norz(ifunc,nn,il,isd,irgm)
227 thetd(i+1)=thetd(i)+dthe
239 xxz1=(yyz1*xlamb)/(2.0d0*
pi)
240 xxz2=(yyz2*xlamb)/(2.0d0*
pi)
241 xxz3=(yyz3*xlamb)/(2.0d0*
pi)
251 x=(2.0d0*
pi*rbar(i))/xlamb
255 call dbmie(x,rfr1,rfi1,thetd,jx,qext1,qscat1,ctbrqs,
257 factor=wt1*
pi*rbar(i)**2*1.0e-8
258 qext=qext+qext1*factor
259 qscat=qscat+qscat1*factor
266 if(rbar(i).ge.rmin1 .and. rbar(i).le.rmax1)
then
267 call dbmie(x,rfr1,rfi1,thetd,jx,qext1,qscat1,ctbrqs,
269 factor=wt1*
pi*rbar(i)**2*1.0e-8
270 qext=qext+qext1*factor
271 qscat=qscat+qscat1*factor
274 if(rbar(i).ge.rmin2 .and. rbar(i).le.rmax2)
then
275 call dbmie(x,rfr2,rfi2,thetd,jx,qext1,qscat1,ctbrqs,
277 factor=wt2*
pi*rbar(i)**2*1.0e-8
278 qext=qext+qext1*factor
279 qscat=qscat+qscat1*factor
288 if(rbar(i).ge.rmin1 .and. rbar(i).le.rmax1)
then
289 call dbmie(x,rfr1,rfi1,thetd,jx,qext1,qscat1,ctbrqs,
291 factor=wt1*
pi*rbar(i)**2*1.0e-8
292 qext=qext+qext1*factor
293 qscat=qscat+qscat1*factor
296 if(rbar(i).ge.rmin2 .and. rbar(i).le.rmax2)
then
297 call dbmie(x,rfr2,rfi2,thetd,jx,qext1,qscat1,ctbrqs,
299 factor=wt2*
pi*rbar(i)**2*1.0e-8
300 qext=qext+qext1*factor
301 qscat=qscat+qscat1*factor
304 if(rbar(i).ge.rmin3 .and. rbar(i).le.rmax3)
then
305 call dbmie(x,rfr3,rfi3,thetd,jx,qext1,qscat1,ctbrqs,
307 factor=wt3*
pi*rbar(i)**2*1.0e-8
308 qext=qext+qext1*factor
309 qscat=qscat+qscat1*factor
316 const=((alambd**2)/(4.0*
pi**2*qscat)) * 4.0*
pi
318 rr(i,1)=t(i,1) *
const
319 rr(i,2)=t(i,2) *
const
320 rr(i,3)=t(i,3) *
const
321 rr(i,4)=t(i,4) *
const
322 if(i .le. jx) xx=thetd(i)
323 if(i .gt. jx) xx=180.0d0-thetd(1801-i+1)
327 angl(i) = dfloat(i-1)/10.0d0
328 if(i .le. jx) xx=(thetd(i)*
pi)/180.0d0
329 if(i .gt. jx) xx=((180.0d0-thetd(1801-i+1))*
pi)/180.0d0
331 rr(i,1)=
const*(t(i,1)+t(i,2)*y**2+2.0d0*t(i,3)*y)
332 rr(i,2)=
const*(t(i,2)+t(i,1)*y**2+2.0d0*t(i,3)*y)
333 rr(i,3)=
const*((t(i,1)+t(i,2))*y+t(i,3)*(1.0d0+y**2))
334 rr(i,4)=
const*(1.d0-y**2)*t(i,4)
335 phfu(i)=(rr(i,1)+rr(i,2))/2.0
336 pol(i)=(rr(i,2)-rr(i,1))/(rr(i,1)+rr(i,2))
346 goto(170,172,174,176,178),ifunc
349 write(8,249)ifunc,mfunc
350 write(8,251)xlamb,rfr1,rfi1
351 write(8,252)ddd(10),ddd(11)
352 write(8,253)xxz1,yyz1
353 write(8,270)ddd(1),ddd(2),ddd(3)
355 write(7,249)ifunc,mfunc
356 write(7,251)xlamb,rfr1,rfi1
357 write(7,252)ddd(10),ddd(11)
358 write(7,253)xxz1,yyz1
359 write(7,270)ddd(1),ddd(2),ddd(3)
363 write(8,249)ifunc,mfunc
364 write(8,251)xlamb,rfr1,rfi1
365 write(8,252)ddd(10),ddd(11)
366 write(8,253)xxz1,yyz1
367 write(8,272)ddd(1),ddd(2),ddd(3),ddd(4)
369 write(7,249)ifunc,mfunc
370 write(7,251)xlamb,rfr1,rfi1
371 write(7,252)ddd(10),ddd(11)
372 write(7,253)xxz1,yyz1
373 write(7,272)ddd(1),ddd(2),ddd(3),ddd(4)
377 write(8,249)ifunc,mfunc
378 write(8,251)xlamb,rfr1,rfi1
379 write(8,252)ddd(10),ddd(11)
380 write(8,253)xxz1,yyz1
382 write(8,275)ddd(1),ddd(2),ddd(3)
384 write(8,274)ddd(1),ddd(2),ddd(3)
387 write(7,249)ifunc,mfunc
388 write(7,251)xlamb,rfr1,rfi1
389 write(7,252)ddd(10),ddd(11)
390 write(7,253)xxz1,yyz1
392 write(7,275)ddd(1),ddd(2),ddd(3)
394 write(7,274)ddd(1),ddd(2),ddd(3)
399 write(8,249)ifunc,mfunc
400 write(8,254)xlamb,rfr1,rfi1,rfr2,rfi2
401 write(8,255)rmin1,rmax1,rmin2,rmax2
402 write(8,256)xxz1,yyz1,xxz2,yyz2
404 write(8,277)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
406 write(8,276)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
409 write(7,249)ifunc,mfunc
410 write(7,254)xlamb,rfr1,rfi1,rfr2,rfi2
411 write(7,255)rmin1,rmax1,rmin2,rmax2
412 write(7,256)xxz1,yyz1,xxz2,yyz2
414 write(7,277)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
416 write(7,276)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
421 write(8,249)ifunc,mfunc
422 write(8,257)xlamb,rfr1,rfi1,rfr2,rfi2,rfr3,rfi3
423 write(8,258)rmin1,rmax1,rmin2,rmax2,rmin3,rmax3
424 write(8,259)xxz1,yyz1,xxz2,yyz2,xxz3,yyz3
426 write(8,279)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
427 1 ddd(7),ddd(8),ddd(9)
429 write(8,278)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
430 1 ddd(7),ddd(8),ddd(9)
433 write(7,249)ifunc,mfunc
434 write(7,257)xlamb,rfr1,rfi1,rfr2,rfi2,rfr3,rfi3
435 write(7,258)rmin1,rmax1,rmin2,rmax2,rmin3,rmax3
436 write(7,259)xxz1,yyz1,xxz2,yyz2,xxz3,yyz3
438 write(7,279)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
439 1 ddd(7),ddd(8),ddd(9)
441 write(7,278)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
442 1 ddd(7),ddd(8),ddd(9)
445 write(8,281)r11,r22,r33,r44,reff,veff
446 write(8,280)ccnsml,bsr,salb,asf,qscat,qext
447 write(7,281)r11,r22,r33,r44,reff,veff
448 write(7,280)ccnsml,bsr,salb,asf,qscat,qext
451 write(8,49)
angl(i),phfu(i),pol(i)
452 write(7,47) (t(i,k),k=1,4),
angl(i)
457 47
format(1p,2d18.12,1p,2d19.12,0p,f6.2)
458 49
format(f5.1,1p,2e15.5)
459 249
format(t2,
'ifunc',t16,
'ifunc_sub'/i5,11x,i5)
460 251
format(t2,
'wavelength',t16,
'refr',t27,
'refi'/1p3e11.3)
461 252
format(t5,
'rmin',t16,
'rmax',t24/1p2e11.3)
462 253
format(t4,
'deltar',t15,
'deltax'/1p2e11.3)
465 254
format(t2,
'wavelength',t16,
'refr1',t27,
'refi1',t38,
'refr2',
466 1 t49,
'refi2'/1p5e11.4)
467 255
format(t5,
'rmin1',t16,
'rmax1',t27,
'rmin2',t38,
'rmax2'/
469 256
format(t4,
'deltar1',t15,
'deltax1',t26,
'deltar2',
470 1 t37,
'deltax2'/1p4e11.3)
471 257
format(t2,
'wavelength',t15,
'refr1',t26,
'refi1',t27,
'refr2',
472 1 t38,
'refi2',t49,
'refr3',t60,
'refi3'/1p7e11.3)
473 258
format(t5,
'rmin1',t16,
'rmax1',t27,
'rmin2',t38,
'rmax2't49,
'rmin3',
474 1 t60,
'rmax3',t68/1p6e11.3)
475 259
format(t4,
'deltar1',t15,
'deltax1',t26,
'deltar2',t37,
'deltax2',
476 1 t48,
'deltar3',t59,
'deltax3'/1p6e11.3)
478 271
format(t2,
'wavelength',t16,
'refr',t27,
'refi't38,
'rmin',
479 1 t49,
'rmax',t51,
'deltar',t62,
'deltax'/1p7e11.3)
480 270
format(t5,
'rc',t16,
'c',t27,
'nu'/1p3e11.3)
481 272
format(t5,
'aa',t16,
'b',t27,
'alfa',t38,
'gama'/1p4e11.3)
482 274
format(t5,
'rg1',t16,
'sig1',t27,
'num1'/1p3e11.3)
483 275
format(t5,
'rg1',t15,
'sig1(ln)',t27,
'num1'/1p3e11.3)
484 276
format(t5,
'rg1't16,
'sig1',t27,
'num1',t38,
'rg2',
485 1 t49,
'sig2',t61,
'num2'/1p6e11.3)
486 277
format(t5,
'rg1't15,
'sig1(ln)',t27,
'num1',t38,
'rg2',
487 1 t48,
'sig2(ln)',t61,
'num2'/1p6e11.3)
488 278
format(t3,
'rg1',t11,
'sig1',t19,
'num1',t27,
'rg2',t35,
'sig2',
489 1 t43,
'num2',t51,
'rg3',t59,
'sig3',t67,
'num3'/9f8.5)
490 279
format(t3,
'rg1',t10,
'sig1(ln)',t19,
'num1',t27,
'rg2',
491 1 t34,
'sig2(ln)',t43,
'num2',t51,
'rg3',t58,
'sig3(ln)',
493 280
format(t4,
'ccnsml',t16,
'bsr',t27,
'salb',t38,
'asf',
494 1 t49,
'qscat',t61,
'qext'/1p4e11.3,1p2e12.4)
495 281
format(t5,
'r11',t16,
'r22',t27,
'r33',t38,
'r44',
496 1 t49,
'reff't60,
'veff'/1p6e11.3)
497 282
format(t3,
'angle',t11,
'phs_func',t26,
'deg_pol')
503 subroutine norz(ifunc,nn,il,isd,irgm)
512 implicit real *8 (a-h,o-z)
514 include
'common_phs.cmn'
525 goto(25,30,35,40,45),ifunc
534 if(rmax.gt.ddd(15))rmax=ddd(15)
540 write(4,710)xlamb,ddd(10),ddd(11)
541 write(4,712)refr1,refi1
542 write(4,720)ddd(1),ddd(2),ddd(3)
553 if(rmax.gt.ddd(15))rmax=ddd(15)
559 write(4,710)xlamb,ddd(10),ddd(11)
560 write(4,712)refr1,refi1
561 write(4,730)ddd(1),ddd(2),ddd(3),ddd(4)
573 rmin1=dexp(dlog(rm1)-4.0d0*sig11)
574 if(rmin1 .lt. 0.00d0)rmin1=0.00d0
575 rmax1=dexp(dlog(rm1)+4.0d0*sig11)
576 if(rmax1.gt.ddd(15))rmax1=ddd(15)
582 write(4,710)xlamb,ddd(10),ddd(11)
583 write(4,712)refr1,refi1
585 write(4,741)ddd(1),ddd(2),ddd(3)
587 write(4,740)ddd(1),ddd(2),ddd(3)
606 rmin1=dexp(dlog(rm1)-4.0d0*sig11)
607 rmax1=dexp(dlog(rm1)+4.0d0*sig11)
608 if(rmin1 .lt. 0.0d0)rmin=0.0d0
609 rmin2=dexp(dlog(rm2)-4.0d0*sig21)
610 rmax2=dexp(dlog(rm2)+4.0d0*sig21)
611 if(rmax1.gt.ddd(15))rmax1=ddd(15)
612 if(rmax2.gt.ddd(15))rmax2=ddd(15)
618 write(4,710)xlamb,ddd(10),ddd(11)
619 write(4,714)refr1,refi1,refr2,refi2
621 write(4,751)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
623 write(4,750)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6)
648 rmin1=dexp(dlog(rm1)-4.0d0*sig11)
649 rmax1=dexp(dlog(rm1)+4.0d0*sig11)
650 if(rmin1 .lt. 0.0d0)rmin=0.0d0
651 rmin2=dexp(dlog(rm2)-4.0d0*sig21)
652 rmax2=dexp(dlog(rm2)+4.0d0*sig21)
653 rmin3=dexp(dlog(rm3)-4.0d0*sig31)
654 rmax3=dexp(dlog(rm3)+4.0d0*sig31)
655 if(rmax1.gt.ddd(15))rmax1=ddd(15)
656 if(rmax2.gt.ddd(15))rmax2=ddd(15)
657 if(rmax3.gt.ddd(15))rmax3=ddd(15)
663 write(4,710)xlamb,ddd(10),ddd(11)
664 write(4,716)refr1,refi1,refr2,refi2,refr3,refi3
666 write(4,760)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
667 1 ddd(7),ddd(8),ddd(9)
669 write(4,760)ddd(1),ddd(2),ddd(3),ddd(4),ddd(5),ddd(6),
670 1 ddd(7),ddd(8),ddd(9)
675 if(ifunc .gt. 1)
goto 53
676 delr=(xlamb*xdx(isd,il,1))/(2.0d0*
pi)
677 n1=(
rc-rmin)/delr+1.0d0
678 n2=(rmax-
rc)/delr+1.0d0
679 xx1=(
rc-rmin)/dfloat(n1)
680 xx2=(rmax-
rc)/dfloat(n2)
686 r(i+1+n1)=r(i+n1)+xx2
694 delr=(xlamb*xdx(isd,il,1))/(2.0d0*
pi)
695 nn=(rmax-rmin)/delr+1.0d0
696 xxx=(rmax-rmin)/dfloat(nn)
700 delr1=(xlamb*xdx(isd,il,1))/(2.0d0*
pi)
701 delr2=(xlamb*xdx(isd,il,2))/(2.0d0*
pi)
702 nn1=(rmax1-rmin1)/delr1+1.0d0
703 xxx1=(rmax1-rmin1)/dfloat(nn1)
704 nn2=(rmax2-rmax1)/delr2+1.0d0
705 xxx2=(rmax2-rmax1)/dfloat(nn2)
709 delr1=(xlamb*xdx(isd,il,1))/(2.0d0*
pi)
710 delr2=(xlamb*xdx(isd,il,2))/(2.0d0*
pi)
711 delr3=(xlamb*xdx(isd,il,3))/(2.0d0*
pi)
712 nn1=(rmax1-rmin1)/delr1+1.0d0
713 xxx1=(rmax1-rmin1)/dfloat(nn1)
714 nn2=(rmax2-rmax1)/delr2+1.0d0
715 xxx2=(rmax2-rmax1)/dfloat(nn2)
716 nn3=(rmax3-rmax2)/delr3+1.0d0
717 xxx3=(rmax3-rmax2)/dfloat(nn3)
725 if(r(i+1).gt.ddd(15))
goto 56
733 if(r(i+1).gt.ddd(15))
goto 56
742 if(r(i+1).gt.ddd(15))
goto 56
750 goto(60,65,70,75,77),ifunc
755 rbar(i)=(r(i+1)+r(i) )/2.0d0
757 if(r(i+1) .gt.
rc)
goto 62
761 sumn1(i+1)=sumn1(i)+wt(i)
764 dn1(i)=c*(
rc/rbar(i))**(xnu+1)
765 dnz1(i)=c*rbar(i)*(
rc/rbar(i))**(xnu+1)
767 sumn1(i+1)=sumn1(i)+wt(i)
769 write(4,118)i,rbar(i),dnz1(i),dn1(i),wt(i),sumn1(i)
776 rbar(i)=(r(i+1)+r(i) )/2.0d0
778 dn1(i)=a*rbar(i)**alfa*dexp(-b*rbar(i)**gama )
779 dnz1(i)=rbar(i)*dn1(i)
781 sumn1(i+1)=sumn1(i)+wt(i)
782 write(4,118)i,rbar(i),dnz1(i),dn1(i),wt(i),sumn1(i)
789 rbar(i)=(r(i+1)+r(i) )/2.0d0
799 consr1=par1*1.0d0/(dsqrt(2.0d0*
pi)*xsig11)
800 dnz1(i)=consr1*dexp(-0.5d0*((xlgr-xlgrm1)/xsig11)**2)
801 dv1(i)=(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz1(i)
802 dn1(i)=(dnz1(i)/rbar(i))
804 sumn1(i+1)=sumn1(i)+wt(i)
807 write(4,118)i,rbar(i),dnz1(i),dn1(i),dv1(i),wt(i),sumn1(i)
818 rbar(i)=(r(i+1)+r(i) )/2.0d0
830 consr1=1.0d0/(dsqrt(2.0d0*
pi)*xsig11)
831 consr2=1.0d0/(dsqrt(2.0d0*
pi)*xsig21)
832 dnz1(i)=consr1*dexp(-0.5d0*((xlgr-xlgrm1)/xsig11)**2)
833 dnz2(i)=consr2*dexp(-0.5d0*((xlgr-xlgrm2)/xsig21)**2)
834 dn1(i)=dnz1(i)/rbar(i)
835 dn2(i)=dnz2(i)/rbar(i)
836 dv1(i)=(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz1(i)
837 dv2(i)=(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz2(i)
843 dnzp(i)=dnz1(i)+dnz2(i)
844 dvp(i)=par1*dv1(i)+par2*dv2(i)
845 wt(i)=dn1(i)*xxx+dn2(i)*xxx
846 sumn1(i+1)=sumn1(i)+dn1(i)*xxx
847 sumn2(i+1)=sumn2(i)+dn2(i)*xxx
848 sumnp(i+1)=sumn1(i+1)+sumn2(i+1)
849 write(4,121)i,rbar(i),dnzp(i),dnp(i),dvp(i),wt(i),sumnp(i)
862 rbar(i)=(r(i+1)+r(i) )/2.0d0
877 consr1=par1/(dsqrt(2.0d0*
pi)*xsig11)
878 consr2=par2/(dsqrt(2.0d0*
pi)*xsig21)
879 consr3=par3/(dsqrt(2.0d0*
pi)*xsig31)
880 dnz1(i)=consr1*dexp(-0.5d0*((xlgr-xlgrm1)/xsig11)**2)
881 dnz2(i)=consr2*dexp(-0.5d0*((xlgr-xlgrm2)/xsig21)**2)
882 dnz3(i)=consr3*dexp(-0.5d0*((xlgr-xlgrm3)/xsig31)**2)
883 dnzp(i)=dnz1(i)+dnz2(i)+dnz3(i)
884 dv1(i)=par1*(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz1(i)
885 dv2(i)=par2*(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz2(i)
886 dv3(i)=par3*(4.0d0/3.0d0)*
pi*rbar(i)**3*dnz3(i)
887 dvp(i)=dv1(i)+dv2(i)+dv3(i)
888 dn1(i)=(dnz1(i)/rbar(i))
889 dn2(i)=(dnz2(i)/rbar(i))
890 dn3(i)=(dnz3(i)/rbar(i))
891 dnp(i)=dn1(i)+dn2(i)+dn3(i)
892 wt(i)=dn1(i)*xxx+dn2(i)*xxx+dn3(i)*xxx
893 sumn1(i+1)=sumn1(i)+dn1(i)*xxx
894 sumn2(i+1)=sumn2(i)+dn2(i)*xxx
895 sumn3(i+1)=sumn3(i)+dn3(i)*xxx
896 sumnp(i+1)=sumn1(i+1)+sumn2(i+1)+sumn3(i+1)
897 write(4,123)i,rbar(i),dnzp(i),dnp(i),dvp(i),wt(i),sumnp(i)
902 118
format(i4,1x,1p6e11.3)
903 121
format(i4,1x,1p6e11.3)
904 123
format(i4,1x,1p6e11.3)
906 710
format(t4,
'lambda',t16,
'rmin',t27,
'rmax'/1p3e11.3)
907 712
format(t5,
'refr1',t16,
'refi1'/1p2e11.3)
908 714
format(t5,
'refr1',t16,
'refi1',t27,
'refr2',t38,
'refi2'/
909 1 1p2e11.3,1x,1p2e11.3)
910 716
format(t5,
'refr1',t16,
'refi1',t27,
'refr2',t38,
'refi2',t49,
'refr3',
911 1 t60,
'refi3'/1p2e11.3,1x,1p2e11.3,1x,1p2e11.3)
912 720
format(t5,
'rc',t16,
'c',t27,
'nu_star'/1p3e11.3)
913 730
format(t5,
'aa',t16,
'b',t27,
'alfa',t38,
'gama'/1p4e11.3)
914 740
format(t5,
'rg1',t16,
'sig1',t27,
'num1'/1p3e11.3)
915 741
format(t5,
'rg1',t15,
'sig1(ln)',t27,
'num1'/1p3e11.3)
916 750
format(t5,
'rg1',t16,
'sig1',t27,
'num1',t38,
'rg2',t49,
'sig2',t60,
918 751
format(t5,
'rg1',t15,
'sig1(ln)',t27,
'num1',t38,
'rg2',t48,
919 1
'sig2(ln)',t60,
'num2'/1p6e11.3)
920 760
format(t4,
'rg1',t12,
'sig1',t20,
'num1',t28,
'rg2',t36,
'sig2',
921 1 t44,
'num2',t52,
'rg3',t60,
'sig3',t68,
'num3'/9f8.5)
922 761
format(t4,
'rg1',t11,
'sig1(ln)',t20,
'num1',t28,
'rg2',t35,
923 1
'sig2(ln)'t44,
'num2',t52,
'rg3',t59,
'sig3(ln)',
925 800
format(t2,
'num',t10,
'rbar',t21,
'dnz',t32,
'dn',t43,
'wt',
927 810
format(t3,
'num',t10,
'rbar',t19,
'dn/dlogr',t32,
'dn/dr',t43,
'wt',
929 820
format(t3,
'num',t10,
'rbar',t19,
'dn/dlogr',t32,
'dn/dr',
930 1 t41,
'dv/dlogr',t54,
'wt',t64,
'tot_num')
931 830
format(t3,
'num',t10,
'rbar',t19,
'dn/dlogr',t32,
'dn/dr',
932 1 t41,
'dv/dlogr',t54,
'wt',t64,
'tot_num')
933 840
format(t3,
'num',t10,
'rbar',t19,
'dn/dlogr',t32,
'dn/dr',
934 1 t41,
'dv/dlogr',t54,
'wt',t64,
'tot_num')
939 subroutine ccn(ifunc,irgm)
945 implicit real*8 (a-h,o-z)
946 include
'common_phs.cmn'
958 aa=(dlog(r0/rg))/(sigln*dsqrt(2.0d0))
960 ccnsml=0.5d0*(1.0d0-erfa1)
967 subroutine dbmie (x,rfr,rfi,thetd,jx,qext,qscat,ctbrqs,eltrmx,
999 implicit real * 8 (a-h,o-z)
1000 real*8 x,rx,rfr,rfi,qext,qscat,t(5),ta(4),tb(2),tc(2)
1001 real*8 td(2),te(2),ctbrqs
1003 real*8 eltrmx(4,1000,2),
pi(3,1000),tau(3,1000),cstht(1000),
1004 1 si2tht(1000),thetd(1000)
1005 complex*16 rf,rrf,rrfx,wm1,fna,fnb,tc1,tc2,wfn(2),acap(7000)
1006 complex*16 fnap,fnbp
1014 real * 8 tt(1801,4),den,az,bz,rt1,rt2,rt3,rt4
1016 equivalence(fnap, td(1)), (fnbp, te(1))
1018 if ( jx .lt. 1001 )
go to 20
1022 20 rf=dcmplx(rfr,-rfi)
1026 t(1)=(x**2)*(rfr**2+rfi**2)
1028 nmx1 = 1.10d0 * t(1)
1029 if (nmx1 .lt. 7000)
go to 21
1033 if (nmx1 .gt. 150)
go to 22
1036 22 acap(nmx1 + 1 ) = ( 0.0d0, 0.0d0 )
1039 acap(nn) = (nn+1) * rrfx - 1.0d0/((nn+1)*rrfx + acap(nn+1))
1042 if ( thetd(j) .lt. 0.0d0 ) thetd(j) = dabs(thetd(j))
1043 if ( thetd(j) .gt. 0.0d0 )
go to 24
1047 24
if ( thetd(j) .ge. 90.0d0 )
go to 25
1048 t(1) = ( 3.1415926535897932d0*thetd(j))/180.0d0
1049 cstht(j) = dcos(t(1))
1050 si2tht(j) = 1.0d0 - cstht(j)**2
1052 25
if ( thetd(j) .gt. 90.0d0 )
go to 28
1056 28
write (6, 5) thetd(j)
1066 300 tau(2, j) = cstht(j)
1070 wm1=dcmplx(t(1),-t(2))
1071 wfn(1)=dcmplx(t(2),t(1))
1072 wfn(2) = rx * wfn(1) - wm1
1073 tc1 = acap(1) * rrf + rx
1074 tc2 = acap(1) * rf + rx
1075 fna = (tc1 * ta(3) - ta(1)) / (tc1 * wfn(2) - wfn(1))
1076 fnb = (tc2 * ta(3) - ta(1)) / (tc2 * wfn(2) - wfn(1))
1080 tb(1) = t(1) * tb(1)
1081 tb(2) = t(1) * tb(2)
1082 tc(1) = t(1) * tc(1)
1083 tc(2) = t(1) * tc(2)
1085 45
if(cstht(j).eq.1.)goto9798
1086 eltrmx(1,j,1)=tc(1)*
pi(2,j)+tb(1)*tau(2,j)
1087 eltrmx(2,j,1)=tc(2)*
pi(2,j)+tb(2)*tau(2,j)
1088 eltrmx(3,j,1)=tb(1)*
pi(2,j)+tc(1)*tau(2,j)
1089 eltrmx(4,j,1)=tb(2)*
pi(2,j)+tc(2)*tau(2,j)
1090 eltrmx(1,j,2)=tc(1)*
pi(2,j)-tb(1)*tau(2,j)
1091 eltrmx(2,j,2)=tc(2)*
pi(2,j)-tb(2)*tau(2,j)
1092 eltrmx(3,j,2)=tb(1)*
pi(2,j)-tc(1)*tau(2,j)
1093 eltrmx(4,j,2)=tb(2)*
pi(2,j)-tc(2)*tau(2,j)
1095 9798 eltrmx(1,j,1)=pip(2,j)*(tc(1)-tb(1))-
pi(2,j)*tb(1)
1096 eltrmx(2,j,1)=pip(2,j)*(tc(2)-tb(2))-
pi(2,j)*tb(2)
1097 eltrmx(3,j,1)=pip(2,j)*(tb(1)-tc(1))-
pi(2,j)*tc(1)
1098 eltrmx(4,j,1)=pip(2,j)*(tb(2)-tc(2))-
pi(2,j)*tc(2)
1099 eltrmx(1,j,2)=pip(2,j)*(tc(1)+tb(1))-
pi(2,j)*tb(1)
1100 eltrmx(2,j,2)=pip(2,j)*(tc(2)+tb(2))-
pi(2,j)*tb(2)
1101 eltrmx(3,j,2)=pip(2,j)*(tc(1)+tb(1))-
pi(2,j)*tc(1)
1102 eltrmx(4,j,2)=pip(2,j)*(tc(2)+tb(2))-
pi(2,j)*tc(2)
1104 if (j .le. jx)
go to 45
1105 qext = 2.0d0 * ( tb(1) + tc(1))
1106 qscat =(tb(1)**2 + tb(2)**2 + tc(1)**2 + tc(2)**2)/0.75d0
1113 pi(3,j)=(t(1)*
pi(2,j)*cstht(j)-n*
pi(1,j))/t(2)
1114 tau(3,j)=cstht(j)*(
pi(3,j)-
pi(1,j))-t(1)*si2tht(j)*
pi(2,j)+tau(1,j
1116 pip(3,j)=t(1)*
pi(2,j)+pip(1,j)
1120 wfn(2) = t(1) * rx * wfn(1) - wm1
1121 tc1 = acap(n) * rrf + n * rx
1122 tc2 = acap(n) * rf + n * rx
1123 fna = (tc1 *ta(3) - ta(1)) / (tc1 * wfn(2) - wfn(1))
1124 fnb = (tc2 * ta(3) - ta(1)) / (tc2 * wfn(2) - wfn(1))
1126 t(4) = t(1) / (t(5) * t(2))
1127 t(2) = (t(2)*(t(5) + 1.0d0))/t(5)
1128 ctbrqs = ctbrqs + t(2) * (td(1) * tb(1) + td(2) * tb(2) + te(1) *
1129 $tc(1) + te(2) * tc(2)) + t(4) * (td(1) * te(1) + td(2) * te(2))
1130 qext = qext + t(3) * (tb(1) + tc(1))
1131 t(4) = tb(1)**2 + tb(2)**2 + tc(1)**2 + tc(2)**2
1132 qscat = qscat + t(3) * t(4)
1137 if(cstht(j).eq.1.)goto86
1138 eltrmx(1,j,1)=eltrmx(1,j,1)+t(1)*(tc(1)*
pi(3,j)+tb(1)*tau(3,j))
1139 eltrmx(2,j,1)=eltrmx(2,j,1)+t(1)*(tc(2)*
pi(3,j)+tb(2)*tau(3,j))
1140 eltrmx(3,j,1)=eltrmx(3,j,1)+t(1)*(tb(1)*
pi(3,j)+tc(1)*tau(3,j))
1141 eltrmx(4,j,1)=eltrmx(4,j,1)+t(1)*(tb(2)*
pi(3,j)+tc(2)*tau(3,j))
1144 eltrmx(1,j,2)=eltrmx(1,j,2)+t(1)*az*(tc(1)*
pi(3,j)-tb(1)*tau(3,j))
1145 eltrmx(2,j,2)=eltrmx(2,j,2)+t(1)*az*(tc(2)*
pi(3,j)-tb(2)*tau(3,j))
1146 eltrmx(3,j,2)=eltrmx(3,j,2)+t(1)*az*(tb(1)*
pi(3,j)-tc(1)*tau(3,j))
1147 eltrmx(4,j,2)=eltrmx(4,j,2)+t(1)*az*(tb(2)*
pi(3,j)-tc(2)*tau(3,j))
1149 86 eltrmx(1,j,1)=eltrmx(1,j,1)+t(1)*(pip(3,j)*(tc(1)-tb(1))-
pi(3,j)*
1151 eltrmx(2,j,1)=eltrmx(2,j,1)+t(1)*(pip(3,j)*(tc(2)-tb(2))-
pi(3,j)*
1153 eltrmx(3,j,1)=eltrmx(3,j,1)+t(1)*(pip(3,j)*(tb(1)-tc(1))-
pi(3,j)*
1155 eltrmx(4,j,1)=eltrmx(4,j,1)+t(1)*(pip(3,j)*(tb(2)-tc(2))-
pi(3,j)*
1162 79 eltrmx(1,j,2)=eltrmx(1,j,2)+t(1)*(az*pip(3,j)*(tc(1)+tb(1))-bz*
1164 eltrmx(2,j,2)=eltrmx(2,j,2)+t(1)*(az*pip(3,j)*(tc(2)+tb(2))-bz*
1166 eltrmx(3,j,2)=eltrmx(3,j,2)+t(1)*(az*pip(3,j)*(tc(1)+tb(1))-bz*
1168 eltrmx(4,j,2)=eltrmx(4,j,2)+t(1)*(az*pip(3,j)*(tc(2)+tb(2))-bz*
1171 if( t(4) .lt. 1.0d-14 )
go to 100
1176 tau(1, j) = tau(2, j)
1177 tau(2, j) = tau(3, j)
1183 if (n .le. nmx2)
go to 65
1191 115 t(i)=eltrmx(i,j,k)
1192 if(cstht(j).eq.1.)goto125
1195 rt3=(t(3)-az*cstht(j)*t(1))/den
1196 rt4=(t(4)-az*cstht(j)*t(2))/den
1197 rt1=(t(1)-az*cstht(j)*t(3))/den
1198 rt2=(t(2)-az*cstht(j)*t(4))/den
1205 if( k .eq. 2) m=1802-m
1206 if(k .eq. 2 .and. m .eq. 901 )
go to 120
1207 tt(m,1) = tt(m,1)+(rt1**2+rt2**2)*wat
1208 tt(m,2) = tt(m,2)+(rt3**2+rt4**2)*wat
1209 tt(m,3) = tt(m,3)+(rt1*rt3+rt2*rt4)*wat
1210 tt(m,4) = tt(m,4)+(rt1*rt4-rt2*rt3)*wat
1212 t(1) = 2.0d0 * rx**2
1214 qscat = qscat * t(1)
1215 ctbrqs = 2.0d0 * ctbrqs * t(1)
1218 5
format(10x
' the value of the scattering angle is greater than 90.0
1219 $ degrees. it is ',e15.4)
1220 6
format(//10x,
' please read the comments'//)
1221 7
format(//10x,
' the value of the argument jx is greater the 1000')
1222 8
format(//10x,
'the upper limit for acap is not enough. suggest get
1223 1detailed output and modify subroutine'//)
1227 subroutine asym(x,y,f,bsr)
1236 implicit real*8 (a-h,o-z)
1238 real*8 xmu(m),dxmu(m),xy(m),xb(m)
1244 xmu(i)=cos(x(i)*conv)
1249 xb(i)=dacos(xmu(i))*y(i)
1254 dxmu(i)=dabs(xmu(i)-xmu(i+1))
1261 sumn=sumn+0.5d0*(xy(i)+xy(i+1))*dxmu(i)
1262 sumd=sumd+0.5d0*(y(i)+y(i+1))*dxmu(i)
1263 sumb=sumb+0.5d0*(xb(i)+xb(i+1))*dxmu(i)
1269 bsr=(2.0d0*sumb)/xden
1272 100
format(
'asym sumn,sumd,sumb,xnum,xden,f,bsr',1p7e11.3)
1277 subroutine convtc (num,nchar, loc)
1318 integer istart, irem, itemp
1333 do 100 i = nchar, istart, -1
1335 loc(i) = char(irem - itemp * 10 + ichar(
'0'))
1342 if (irem.ne.0) loc(1) =
'*'
1347 subroutine spline(s,x,y,n,in,t,il,iu,vl,vu,e,u)
1355 implicit real*8 (a-h,o-z)
1356 dimension x(*), y(*), e(*),u(*)
1367 u(1)=(c1-vl)/2.0d0/b1
1382 20 u(j)=(d-c*u(j-1))/p
1385 e(n)=u(n1)/(1.0d0-e(n1))
1387 24 c2=vu-(y(n)-y(n1))/(x(n)-x(n1))
1388 e(n)=(c2/(x(n)-x(n1))-u(n1))/(2.0d0+e(n1))
1394 e(k)=e(k)*e(k+1)+u(k)
1399 u(k)=(y(k+1)-y(k))/b2-b2*(e(k)*2+e(k+1))
1402 u(n)=(e(n1)+2.0d0*e(n))*b2+(y(n)-y(n1))/b2
1403 40
if (x(1).gt.x(n))
go to 50
1411 60
if (s.ge.x(mub+idir))
go to 100
1412 if (s.le.x(mlb+1-idir))
go to 110
1416 70
if (iabs(mu-ml).le.1)
go to 120
1418 if (s.lt.x(mav))
go to 90
1425 110 mu=mlb+2*(1-idir)
1429 t=(e(mu-1)*((x(mu)-s)**3)+e(mu)*((s-x(mu-1))**3)+
1430 1 (y(mu-1)-e(mu-1)*((x(mu)-x(mu-1))**2))*(x(mu)-s)+
1431 2 (y(mu)-e(mu)*((x(mu)-x(mu-1))**2))*(s-x(mu-1)))/
1438 implicit real*8 (a-h,o-z)
1448 implicit real*8 (a-h,o-z)
1449 if(x.lt.0.d0.or.a.le.0.d0) stop
1451 call gser(gamser,a,x,gln)
1454 call gcf(gammcf,a,x,gln)
1460 subroutine gser(gamser,a,x,gln)
1461 implicit real*8 (a-h,o-z)
1476 if(dabs(del).lt.dabs(sum)*eps)
go to 1
1478 stop
'a too large, itmax too small'
1479 1 gamser=sum*dexp(-x+a*dlog(x)-gln)
1483 subroutine gcf(gammcf,a,x,gln)
1484 implicit real*8 (a-h,o-z)
1504 if(dabs((g-gold)/g).lt.eps)
go to 1
1508 stop
'a too large, itmax too small'
1510 gammcf=dexp(-x+a*dlog(x)-gln)*g
1515 implicit real*8 (a-h,o-z)
1516 real*8 cof(6),stp,half,one,fpf,x,tmp,ser
1517 data cof,stp/76.18009173d0,-86.50532033d0,24.01409822d0,
1518 * -1.231739516d0,.120858003d-2,-.536382d-5,2.50662827465d0/
1519 data half,one,fpf/0.5d0,1.0d0,5.5d0/
1522 tmp=(x+half)*log(tmp)-tmp
1533 subroutine moment(ifunc,nn)
1537 implicit real *8 (a-h,o-z)
1539 include
'common_phs.cmn'
1557 if(rbar(i).gt.rmin1 .and. rbar(i).le.rmax1)
then
1558 prod1=rbar(i)*dn1(i)*deltar(i)+prod1
1559 prod2=rbar(i)**2*dn1(i)*deltar(i)+prod2
1560 prod3=rbar(i)**3*dn1(i)*deltar(i)+prod3
1561 prod4=rbar(i)**4*dn1(i)*deltar(i)+prod4
1564 if(ifunc.eq.4 .or. ifunc.eq.5)
then
1566 if(rbar(i).gt.rmin2 .and. rbar(i).le.rmax2)
then
1567 prod1=rbar(i)*dn2(i)*deltar(i)+prod1
1568 prod2=rbar(i)**2*dn2(i)*deltar(i)+prod2
1569 prod3=rbar(i)**3*dn2(i)*deltar(i)+prod3
1570 prod4=rbar(i)**4*dn2(i)*deltar(i)+prod4
1576 if(rbar(i).gt.rmin3 .and. rbar(i).le.rmax3)
then
1577 prod1=rbar(i)*dn3(i)*deltar(i)+prod1
1578 prod2=rbar(i)**2*dn3(i)*deltar(i)+prod2
1579 prod3=rbar(i)**3*dn3(i)*deltar(i)+prod3
1580 prod4=rbar(i)**4*dn3(i)*deltar(i)+prod4
1585 veffz=prod4/(reffz**2*prod2)
1593 100
format(
'prod1,prod2,prod3,prod4,reffz,veffz'/1p6e11.3)
1598 subroutine xmntst(ifunc,irgm)
1605 implicit real*8 (a-h,o-z)
1606 include
'common_phs.cmn'
1607 real*8 q11(3),q22(3),q33(3),q44(3),qnn(3)
1608 real*8 amin(3),amax(3),agm(3),asigln(3)
1612 write(*,*)
'welcome to subroutine xmntst'
1623 aa=(dlog(r0/rg))/(sigln*dsqrt(2.0d0))
1625 ccnsml=0.5d0*(1.0d0-erfa1)
1648 asigln(1)=dlog(ddd(2))
1649 asigln(2)=dlog(ddd(5))
1650 asigln(3)=dlog(ddd(8))
1652 if(nmode.gt.3)
return
1655 q11(i)=agm(i)*dexp(0.5d0*asigln(i)**2)
1656 q22(i)=agm(i)**2*dexp(0.5d0*(2.0d0*asigln(i))**2)
1657 q33(i)=agm(i)**3*dexp(0.5d0*(3.0d0*asigln(i))**2)
1658 q44(i)=agm(i)**4*dexp(0.5d0*(4.0d0*asigln(i))**2)
1662 rgp1=dlog(agm(i))+asigln(i)**2
1663 rgp2=dlog(agm(i))+2.0d0*asigln(i)**2
1664 rgp3=dlog(agm(i))+3.0d0*asigln(i)**2
1665 rgp4=dlog(agm(i))+4.0d0*asigln(i)**2
1667 erf1=erf((dlog(amax(i))-rgp1)/(asigln(i)*dsqrt(2.0d0)))
1668 erf2=erf((dlog(amax(i))-rgp2)/(asigln(i)*dsqrt(2.0d0)))
1669 erf3=erf((dlog(amax(i))-rgp3)/(asigln(i)*dsqrt(2.0d0)))
1670 erf4=erf((dlog(amax(i))-rgp4)/(asigln(i)*dsqrt(2.0d0)))
1672 grf1=erf((dlog(amin(i))-rgp1)/(asigln(i)*dsqrt(2.0d0)))
1673 grf2=erf((dlog(amin(i))-rgp2)/(asigln(i)*dsqrt(2.0d0)))
1674 grf3=erf((dlog(amin(i))-rgp3)/(asigln(i)*dsqrt(2.0d0)))
1675 grf4=erf((dlog(amin(i))-rgp4)/(asigln(i)*dsqrt(2.0d0)))
1677 q11(i)=0.5d0*q11(i)*(erf1-grf1)
1678 q22(i)=0.5d0*q22(i)*(erf2-grf2)
1679 q33(i)=0.5d0*q33(i)*(erf3-grf3)
1680 q44(i)=0.5d0*q44(i)*(erf4-grf4)
1688 r11=r11+qnn(i)*q11(i)
1689 r22=r22+qnn(i)*q22(i)
1690 r33=r33+qnn(i)*q33(i)
1691 r44=r44+qnn(i)*q44(i)
1696 veff=r44/(reff**2*r22)
1698 write(*,100)r11,r22,r33,r44,reff,veff
1699 100
format(
'final: r11,r22,r33,r44,reff,veff'/1p6e11.3)