7 #include <boost/math/special_functions.hpp>
21 void dkia_(
int *IFAC,
double *
X,
double *
A,
double *DKI,
double *DKID,
int *IERRO);
22 void dlia_(
int *IFAC,
double *
X,
double *
A,
double *DLI,
double *DLID,
int *IERRO);
27 #define jn(order, x) boost::math::cyl_bessel_j(order, x)
29 #define yn(order, x) boost::math::cyl_neumann(order, x)
31 #define in(order, x) boost::math::cyl_bessel_i(order, x)
33 #define kn(order, x) boost::math::cyl_bessel_k(order, x)
34 #define limu(im_order, x) Rossegger::Limu(im_order, x)
35 #define kimu(im_order, x) Rossegger::Kimu(im_order, x)
55 char zeroesfilename[200];
56 sprintf(zeroesfilename,
"rosseger_zeroes_eps%1.0E_a%2.2f_b%2.2f_L%2.2f.root",
epsilon,
a,
b,
L);
57 TFile *fileptr = TFile::Open(zeroesfilename,
"READ");
72 std::cout <<
"Rossegger object initialized as follows:" << std::endl;
73 std::cout <<
" Inner Radius = " <<
a <<
" cm." << std::endl;
74 std::cout <<
" Outer Radius = " <<
b <<
" cm." << std::endl;
75 std::cout <<
" Half Length = " <<
L <<
" cm." << std::endl;
79 std::cout <<
"CheckZeroes(0.01) returned false, exiting" << std::endl;
87 double previous = (this->*
func)(order, xstart);
89 double value = previous;
94 value = (this->*
func)(order, x);
95 if (value == 0) std::cout <<
"hit it exactly! Go buy a lottery ticket!" << std::endl;
102 double slope = (value - previous) / epsilon;
103 double intercept = value - slope *
x;
104 double x0 = -intercept / slope;
105 if (
verbosity > 1) std::cout <<
" " << x0 <<
"," << std::endl;
106 double n0 = (this->*
func)(order, x - epsilon);
107 double n1 = (this->*
func)(order, x + epsilon);
108 if ((n0 < 0 && n1 < 0) || (n0 > 0 && n1 > 0))
110 printf(
"neighbors on both sides have the same sign. Check your function resolution!\n");
118 std::cout <<
"logic break!\n";
125 std::cout <<
"Now filling the Beta[m][n] Array..." << std::endl;
129 if (
verbosity) std::cout <<
"Filling Beta[" <<
m <<
"][n]..." << std::endl;
150 N2mn[
m][
n] *= (jna_over_jnb * jna_over_jnb - 1.0);
152 if (
verbosity > 1) std::cout <<
"m: " <<
m <<
" n: " <<
n <<
" N2[m][n]: " <<
N2mn[
m][
n];
153 double integral = 0.0;
157 for (
double r = a;
r <
b;
r +=
step)
159 double rmnval =
Rmn(
m,
n,
r);
161 integral += rmnval * rmnval *
r *
step;
163 std::cout <<
" Int: " << integral << std::endl;
169 std::cout <<
"Done." << std::endl;
174 std::cout <<
"Now filling the Mu[n][k] Array..." << std::endl;
185 if (
verbosity) std::cout <<
"Filling Mu[" <<
n <<
"][k]..." << std::endl;
195 printf(
"adjacent values are Rnk[mu-epsilon]=%E\tRnk[mu+epsilon]=%E\n",
197 if (
verbosity > 100)
printf(
"values of argument to limu and kimu are %f and %f\n",
211 double integral = 0.0;
218 integral += rnkval * rnkval /
r *
step;
222 std::cout <<
" Int: " << integral << std::endl;
225 if (
verbosity > 1) std::cout <<
"n: " <<
n <<
" k: " <<
k <<
" N2nk[n][k]: " <<
N2nk[
n][
k];
229 std::cout <<
"Done." << std::endl;
242 if (
abs(result) > epsilon)
244 printf(
"(m=%d,n=%d) Jm(x)Ym(lx)-Jm(lx)Ym(x) = %f for x=b*%f\n",
m,
n, result,
Betamn[
m][
n]);
256 if (
abs(result) > epsilon * 100)
258 printf(
"(n=%d,k=%d) limu(npi*a/L)kimu(npi*b/L)-kimu(npi*a/L)kimu(npi*b/L) = %f (>eps*100) for mu=%f\n",
259 n, k, result,
Munk[
n][k]);
325 dlia_(&IFAC, &X, &A, &DLI, &DERR, &IERRO);
337 dkia_(&IFAC, &X, &A, &DKI, &DERR, &IERRO);
343 double lx =
a * x /
b;
346 return jn(m, x) *
yn(m, lx) -
jn(m, lx) *
yn(m, x);
351 if (
verbosity > 100) std::cout <<
"Determine Rmn(" << m <<
"," << n <<
"," << r <<
") = ";
357 if (r < a || r >
b) error = 1;
360 std::cout <<
"Invalid arguments Rmn(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
371 if (
verbosity > 100) std::cout << R << std::endl;
377 if (
verbosity > 100) std::cout <<
"Determine Rmn(" << m <<
"," << n <<
"," << r <<
") = ";
383 if (r < a || r >
b) error = 1;
386 std::cout <<
"Invalid arguments Rmn(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
397 if (
verbosity > 100) std::cout << R << std::endl;
407 if (r < a || r >
b) error = 1;
410 std::cout <<
"Invalid arguments Rmn1(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
430 if (r < a || r >
b) error = 1;
433 std::cout <<
"Invalid arguments Rmn1(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
442 double BetaN_ = (n + 1) *
pi /
L;
443 R =
kn(m, BetaN_ *
a) *
in(m, BetaN_ * r) -
in(m, BetaN_ * a) *
kn(m, BetaN_ * r);
454 if (r < a || r >
b) error = 1;
457 std::cout <<
"Invalid arguments Rmn2(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
477 if (r < a || r >
b) error = 1;
480 std::cout <<
"Invalid arguments Rmn2(" << m <<
"," << n <<
"," << r <<
")" << std::endl;
489 double BetaN_ = (n + 1) *
pi /
L;
490 R =
kn(m, BetaN_ * b) *
in(m, BetaN_ * r) -
in(m, BetaN_ * b) *
kn(m, BetaN_ * r);
501 if (ref < a || ref >
b) error = 1;
502 if (r < a || r > b) error = 1;
505 std::cout <<
"Invalid arguments RPrime(" << m <<
"," << n <<
"," << ref <<
"," << r <<
")" << std::endl;
519 double term1 =
kn(m, BetaN_ * ref) * (
in(m - 1, BetaN_ * r) +
in(m + 1, BetaN_ * r));
520 double term2 =
in(m, BetaN_ * ref) * (
kn(m - 1, BetaN_ * r) +
kn(m + 1, BetaN_ * r));
521 R = BetaN_ / 2.0 * (term1 + term2);
532 if (ref < a || ref >
b) error = 1;
533 if (r < a || r > b) error = 1;
536 std::cout <<
"Invalid arguments RPrime(" << m <<
"," << n <<
"," << ref <<
"," << r <<
")" << std::endl;
548 double BetaN_ = (n + 1) *
pi /
L;
549 double term1 =
kn(m, BetaN_ * ref) * (
in(m - 1, BetaN_ * r) +
in(m + 1, BetaN_ * r));
550 double term2 =
in(m, BetaN_ * ref) * (
kn(m - 1, BetaN_ * r) +
kn(m + 1, BetaN_ * r));
551 R = BetaN_ / 2.0 * (term1 + term2);
559 if (
verbosity > 10)
printf(
"Rnk_for_zeroes called with n=%d,mu=%f\n", n, mu);
565 return limu(mu, betana) *
kimu(mu, betanb) -
kimu(mu, betana) *
limu(mu, betanb);
571 if (
verbosity > 10)
printf(
"Rnk_for_zeroes called with n=%d,mu=%f\n", n, mu);
572 double BetaN_ = (n + 1) *
pi /
L;
573 double betana = BetaN_ *
a;
574 double betanb = BetaN_ *
b;
578 return limu(mu, betana) *
kimu(mu, betanb) -
kimu(mu, betana) *
limu(mu, betanb);
587 if (r < a || r >
b) error = 1;
590 std::cout <<
"Invalid arguments Rnk(" << n <<
"," << k <<
"," << r <<
")" << std::endl;
606 if (r < a || r >
b) error = 1;
609 std::cout <<
"Invalid arguments Rnk(" << n <<
"," << k <<
"," << r <<
")" << std::endl;
613 double BetaN_ = (n + 1) *
pi /
L;
624 if (r < a || r >
b) error = 1;
625 if (phi < 0 || phi > 2 *
pi) error = 1;
626 if (z < 0 || z >
L) error = 1;
627 if (r1 < a || r1 > b) error = 1;
628 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
629 if (z1 < 0 || z1 > L) error = 1;
632 std::cout <<
"Invalid arguments Ez(";
633 std::cout << r <<
",";
634 std::cout << phi <<
",";
635 std::cout << z <<
",";
636 std::cout << r1 <<
",";
637 std::cout << phi1 <<
",";
639 std::cout <<
")" << std::endl;
647 if (
verbosity > 10) std::cout << std::endl
653 if (
verbosity > 10) std::cout <<
" " << term;
654 term *= (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));
655 if (
verbosity > 10) std::cout <<
" " << term;
657 if (
verbosity > 10) std::cout <<
" " << term;
666 if (
verbosity > 10) std::cout <<
" " << term;
668 if (
verbosity > 10) std::cout <<
" " << term <<
" " << G << std::endl;
673 if (
verbosity) std::cout <<
"Ez = " << G << std::endl;
683 if (r < a || r >
b) error = 1;
684 if (phi < 0 || phi > 2 *
pi) error = 1;
685 if (z < 0 || z >
L) error = 1;
686 if (r1 < a || r1 > b) error = 1;
687 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
688 if (z1 < 0 || z1 > L) error = 1;
691 std::cout <<
"Invalid arguments Ez(";
692 std::cout << r <<
",";
693 std::cout << phi <<
",";
694 std::cout << z <<
",";
695 std::cout << r1 <<
",";
696 std::cout << phi1 <<
",";
698 std::cout <<
")" << std::endl;
706 if (
verbosity > 10) std::cout << std::endl
711 double term = 1 / (2.0 *
pi);
712 if (
verbosity > 10) std::cout <<
" " << term;
713 term *= (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));
714 if (
verbosity > 10) std::cout <<
" " << term;
716 if (
verbosity > 10) std::cout <<
" " << term;
719 term *= cosh(
Betamn[m][n] * z) * sinh(
Betamn[m][n] * (L - z1)) / sinh(
Betamn[m][n] * L);
723 term *= -cosh(
Betamn[m][n] * (L - z)) * sinh(
Betamn[m][n] * z1) / sinh(
Betamn[m][n] * L);
726 if (
verbosity > 10) std::cout <<
" " << term;
728 if (
verbosity > 10) std::cout <<
" " << term <<
" " << G << std::endl;
731 if (
verbosity) std::cout <<
"Ez = " << G << std::endl;
745 if (r < a || r >
b) error = 1;
746 if (phi < 0 || phi > 2 *
pi) error = 1;
747 if (z < 0 || z >
L) error = 1;
748 if (r1 < a || r1 > b) error = 1;
749 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
750 if (z1 < 0 || z1 > L) error = 1;
753 std::cout <<
"Invalid arguments Er(";
754 std::cout << r <<
",";
755 std::cout << phi <<
",";
756 std::cout << z <<
",";
757 std::cout << r1 <<
",";
758 std::cout << phi1 <<
",";
760 std::cout <<
")" << std::endl;
772 part = (2 - ((
m == 0) ? 1 : 0)) * cos(
m * (phi - phi1));
773 if (
verbosity > 10)
printf(
"(2 - ((m==0)?1:0))*cos(m*(phi-phi1)); = %f\n", part);
776 if (
verbosity > 10)
printf(
"sin(BetaN[n]*z)*sin(BetaN[n]*z1); = %f\n", part);
793 if (
verbosity) std::cout <<
"Er = " << G << std::endl;
806 if (r < a || r >
b) error = 1;
807 if (phi < 0 || phi > 2 *
pi) error = 1;
808 if (z < 0 || z >
L) error = 1;
809 if (r1 < a || r1 > b) error = 1;
810 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
811 if (z1 < 0 || z1 > L) error = 1;
814 std::cout <<
"Invalid arguments Er(";
815 std::cout << r <<
",";
816 std::cout << phi <<
",";
817 std::cout << z <<
",";
818 std::cout << r1 <<
",";
819 std::cout << phi1 <<
",";
821 std::cout <<
")" << std::endl;
831 double term = 1 / (L *
pi);
832 term *= (2 - ((
m == 0) ? 1 : 0)) * cos(
m * (phi - phi1));
833 double BetaN_ = (
n + 1) * pi / L;
834 term *= sin(BetaN_ * z) * sin(BetaN_ * z1);
844 term /= TMath::BesselI(
m, BetaN_ *
a) * TMath::BesselK(
m, BetaN_ * b) - TMath::BesselI(
m, BetaN_ * b) * TMath::BesselK(
m, BetaN_ * a);
850 if (
verbosity) std::cout <<
"Er = " << G << std::endl;
861 if (r < a || r >
b) error = 1;
862 if (phi < 0 || phi > 2 *
pi) error = 1;
863 if (z < 0 || z >
L) error = 1;
864 if (r1 < a || r1 > b) error = 1;
865 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
866 if (z1 < 0 || z1 > L) error = 1;
869 std::cout <<
"Invalid arguments Ephi(";
870 std::cout << r <<
",";
871 std::cout << phi <<
",";
872 std::cout << z <<
",";
873 std::cout << r1 <<
",";
874 std::cout << phi1 <<
",";
876 std::cout <<
")" << std::endl;
894 term *= -sinh(
Munk[
n][
k] * (pi - (phi - phi1)));
898 term *= sinh(
Munk[
n][
k] * (pi - (phi1 - phi)));
906 if (
verbosity) std::cout <<
"Ephi = " << G << std::endl;
916 if (r < a || r >
b) error = 1;
917 if (phi < 0 || phi > 2 *
pi) error = 1;
918 if (z < 0 || z >
L) error = 1;
919 if (r1 < a || r1 > b) error = 1;
920 if (phi1 < 0 || phi1 > 2 * pi) error = 1;
921 if (z1 < 0 || z1 > L) error = 1;
924 std::cout <<
"Invalid arguments Ephi(";
925 std::cout << r <<
",";
926 std::cout << phi <<
",";
927 std::cout << z <<
",";
928 std::cout << r1 <<
",";
929 std::cout << phi1 <<
",";
931 std::cout <<
")" << std::endl;
945 double BetaN_ = (n + 1) * pi / L;
946 double term = 1 / (L *
r);
947 if (
verbosity) std::cout <<
" 1/L=" << term;
948 term *= sin(BetaN_ * z) * sin(BetaN_ * z1);
949 if (
verbosity) std::cout <<
" *sinsin=" << term;
950 term *=
Rnk_(n, k, r) *
Rnk_(n, k, r1);
951 if (
verbosity) std::cout <<
" *rnkrnk=" << term;
953 if (
verbosity) std::cout <<
" */nnknnk=" << term;
958 term *= -sinh(
Munk[n][k] * (pi - (phi - phi1)));
964 term *= sinh(
Munk[n][k] * (pi - (phi1 - phi)));
968 if (
verbosity) std::cout <<
" *sinh(mu*pi-phi-phi)=" << term;
969 term *= 1 / (sinh(pi *
Munk[n][k]));
973 if (
verbosity) std::cout <<
" /sinh=" << term <<
" G=" << G << std::endl;
976 if (
verbosity) std::cout <<
"Ephi = " << G << std::endl;
984 TFile *output = TFile::Open(destfile,
"RECREATE");
987 TTree *tInfo =
new TTree(
"info",
"Mu[n][k] values");
989 tInfo->Branch(
"order", &ord);
990 tInfo->Branch(
"epsilon", &epsilon);
996 TTree *tmunk =
new TTree(
"munk",
"Mu[n][k] values");
997 tmunk->Branch(
"n", &n);
998 tmunk->Branch(
"k", &k);
999 tmunk->Branch(
"munk", &munk);
1000 tmunk->Branch(
"n2nk", &n2nk);
1001 for (n = 0; n < ord; n++)
1003 for (k = 0; k < ord; k++)
1011 TTree *tbetamn =
new TTree(
"betamn",
"Beta[m][n] values");
1012 tbetamn->Branch(
"m", &m);
1013 tbetamn->Branch(
"n", &n);
1014 tbetamn->Branch(
"betamn", &betamn);
1015 tbetamn->Branch(
"n2mn", &n2mn);
1016 for (m = 0; m < ord; m++)
1018 for (n = 0; n < ord; n++)
1036 TFile *
f = TFile::Open(destfile,
"READ");
1037 printf(
"reading rossegger zeroes from %s\n", destfile);
1038 TTree *tInfo = (TTree *) (f->Get(
"info"));
1040 tInfo->SetBranchAddress(
"order", &ord);
1041 tInfo->SetBranchAddress(
"epsilon", &epsilon);
1043 printf(
"order=%d,epsilon=%f\n", ord, epsilon);
1046 double munk, betamn;
1048 TTree *tmunk = (TTree *) (f->Get(
"munk"));
1049 tmunk->SetBranchAddress(
"n", &n);
1050 tmunk->SetBranchAddress(
"k", &k);
1051 tmunk->SetBranchAddress(
"munk", &munk);
1052 tmunk->SetBranchAddress(
"n2nk", &n2nk);
1053 for (
int i = 0; i < tmunk->GetEntries(); i++)
1060 TTree *tbetamn = (TTree *) (f->Get(
"betamn"));
1061 tbetamn->SetBranchAddress(
"m", &m);
1062 tbetamn->SetBranchAddress(
"n", &n);
1063 tbetamn->SetBranchAddress(
"betamn", &betamn);
1064 tbetamn->SetBranchAddress(
"n2mn", &n2mn);
1065 for (
int i = 0; i < tbetamn->GetEntries(); i++)
1067 tbetamn->GetEntry(i);