Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-18 09:17:53

0001 #include "Rossegger.h"
0002 
0003 #include <TFile.h>
0004 #include <TMath.h>
0005 #include <TTree.h>
0006 
0007 #include <boost/math/special_functions.hpp>  //covers all the special functions.
0008 
0009 #include <algorithm>  // for max
0010 #include <cmath>
0011 #include <cstdlib>  // for exit, abs
0012 #include <format>
0013 #include <fstream>
0014 #include <iostream>
0015 #include <sstream>
0016 #include <string>
0017 
0018 // the limu and kimu terms, that i need to think about a little while longer...
0019 extern "C"
0020 {
0021   void dkia_(int *IFAC, double *X, double *A, double *DKI, double *DKID, int *IERRO);
0022   void dlia_(int *IFAC, double *X, double *A, double *DLI, double *DLID, int *IERRO);
0023 }
0024 //
0025 
0026 // Bessel Function J_n(x):
0027 #define jn(order, x) boost::math::cyl_bessel_j(order, x)
0028 // Bessel (Neumann) Function Y_n(x):
0029 #define yn(order, x) boost::math::cyl_neumann(order, x)
0030 // Modified Bessel Function of the first kind I_n(x):
0031 #define in(order, x) boost::math::cyl_bessel_i(order, x)
0032 // Modified Bessel Function of the second kind K_n(x):
0033 #define kn(order, x) boost::math::cyl_bessel_k(order, x)
0034 #define limu(im_order, x) Rossegger::Limu(im_order, x)
0035 #define kimu(im_order, x) Rossegger::Kimu(im_order, x)
0036 
0037 /*
0038   This is a modified/renamed copy of Carlos and Tom's "Spacecharge" class, modified to use boost instead of fortran routines, and with phi terms added.
0039  */
0040 
0041 Rossegger::Rossegger(double InnerRadius, double OuterRadius, double Rdo_Z, double precision)
0042   : a(InnerRadius)
0043   , b(OuterRadius)
0044   , L(Rdo_Z)
0045   , epsilon(precision)
0046 {
0047   PrecalcFreeConstants();
0048 
0049   // load the greens functions:
0050   std::string zeroesfilename = std::format("rosseger_zeroes_eps{:.0E}_a{:.2f}_b{:.2f}_L{:.2f}.root", epsilon, a, b, L);
0051 
0052   TFile *fileptr = TFile::Open(zeroesfilename.c_str(), "READ");
0053   if (!fileptr)
0054   {  // generate the lookuptable
0055     FindMunk(epsilon);
0056     FindBetamn(epsilon);
0057     SaveZeroes(zeroesfilename);
0058   }
0059   else
0060   {  // load it from a file
0061     fileptr->Close();
0062     LoadZeroes(zeroesfilename);
0063   }
0064 
0065   PrecalcDerivedConstants();
0066 
0067   std::cout << "Rossegger object initialized as follows:" << std::endl;
0068   std::cout << "  Inner Radius = " << a << " cm." << std::endl;
0069   std::cout << "  Outer Radius = " << b << " cm." << std::endl;
0070   std::cout << "  Half  Length = " << L << " cm." << std::endl;
0071 
0072   if (!CheckZeroes(0.01))
0073   {
0074     std::cout << "CheckZeroes(0.01) returned false, exiting" << std::endl;
0075     exit(1);
0076   }
0077   return;
0078 }
0079 
0080 double Rossegger::FindNextZero(double xstart, double localepsilon, int order, double (Rossegger::*func)(int, double))
0081 {
0082   double previous = (this->*func)(order, xstart);
0083   double x = xstart + localepsilon;
0084   double value = previous;
0085 
0086   while (!((value == 0) || (value < 0 && previous > 0) || (value > 0 && previous < 0)))
0087   {
0088     //  Rossegger equation 5.12
0089     value = (this->*func)(order, x);
0090     if (value == 0)
0091     {
0092       std::cout << "hit it exactly!  Go buy a lottery ticket!" << std::endl;
0093     }
0094     if ((value == 0) || (value < 0 && previous > 0) || (value > 0 && previous < 0))
0095     {
0096       // when we go from one sign to the other, we have bracketed the zero
0097       // the following is mathematically equivalent to finding the delta
0098       // between the x of our previous value and the point where we cross the axis
0099       // then returning x0=x_old+delta
0100       double slope = (value - previous) / localepsilon;
0101       double intercept = value - slope * x;
0102       double x0 = -intercept / slope;
0103       if (verbosity > 1)
0104       {
0105         std::cout << " " << x0 << "," << std::endl;
0106       }
0107       double n0 = (this->*func)(order, x - localepsilon);
0108       double n1 = (this->*func)(order, x + localepsilon);
0109       if ((n0 < 0 && n1 < 0) || (n0 > 0 && n1 > 0))
0110       {
0111         std::cout << "neighbors on both sides have the same sign.  Check your function resolution!" << std::endl;
0112       }
0113 
0114       return x0;
0115     }
0116     previous = value;
0117     x += localepsilon;
0118   }
0119   std::cout << "logic break!" << std::endl;
0120   exit(1);
0121 }
0122 
0123 void Rossegger::FindBetamn(double localepsilon)
0124 {
0125   std::cout << "Now filling the Beta[m][n] Array..." << std::endl;
0126   if (verbosity > 5)
0127   {
0128     std::cout << "numberOfOrders= " << NumberOfOrders << std::endl;
0129   }
0130   for (int m = 0; m < NumberOfOrders; m++)
0131   {
0132     if (verbosity)
0133     {
0134       std::cout << "Filling Beta[" << m << "][n]..." << std::endl;
0135     }
0136 
0137     double x = localepsilon;
0138     for (int n = 0; n < NumberOfOrders; n++)
0139     {  //  !!!  Off by one from Rossegger convention  !!!
0140       x = FindNextZero(x, localepsilon, m, &Rossegger::Rmn_for_zeroes);
0141       Betamn[m][n] = x / b;
0142       x += localepsilon;
0143     }
0144   }
0145 
0146   //  Now fill in the N2mn array...
0147   for (int m = 0; m < NumberOfOrders; m++)
0148   {
0149     for (int n = 0; n < NumberOfOrders; n++)
0150     {
0151       //  Rossegger Equation 5.17
0152       //  N^2_mn = 2/(pi * beta)^2 [ {Jm(beta a)/Jm(beta b)}^2 - 1 ]
0153       N2mn[m][n] = 2 / (pi * pi * Betamn[m][n] * Betamn[m][n]);
0154       // N2mn[m][n] *= (jn(m,Betamn[m][n]*a)*jn(m,Betamn[m][n]*a)/jn(m,Betamn[m][n]*b)/jn(m,Betamn[m][n]*b) - 1.0);
0155       double jna_over_jnb = jn(m, Betamn[m][n] * a) / jn(m, Betamn[m][n] * b);
0156       N2mn[m][n] *= (jna_over_jnb * jna_over_jnb - 1.0);
0157       // rcc note!  in eq 5.17, N2nm is set with betamn[m][n], but from context that looks to be a typo.  The order is mn everywhere else
0158       if (verbosity > 1)
0159       {
0160         std::cout << "m: " << m << " n: " << n << " N2[m][n]: " << N2mn[m][n];
0161       }
0162       double step = 0.01;
0163       if (verbosity > 1)
0164       {
0165         double integral = 0.0;
0166         for (double r = a; r < b; r += step)
0167         {
0168           double rmnval = Rmn(m, n, r);
0169 
0170           integral += rmnval * rmnval * r * step;  // Rmn(m,n,r)*Rmn(m,n,r)*r*step;
0171         }
0172         std::cout << " Int: " << integral << std::endl;
0173       }
0174       // N2mn[m][n] = integral;
0175     }
0176   }
0177 
0178   std::cout << "Done." << std::endl;
0179 }
0180 
0181 void Rossegger::FindMunk(double localepsilon)
0182 {
0183   std::cout << "Now filling the Mu[n][k] Array..." << std::endl;
0184   if (verbosity > 5)
0185   {
0186     std::cout << "numberOfOrders= " << NumberOfOrders << std::endl;
0187   }
0188   // We're looking for the zeroes of Rossegger eqn. 5.46:
0189   // R_nk(mu_nk;a,b)=Limu(Beta_n*a)Kimu(Beta_n*b)-Kimu(Beta_n*a)Limu(Beta_n*b)=0
0190   // since a and b are fixed, R_nk is a function solely of mu_nk and n.
0191   // for each 'n' we wish to find the a set of k mu_n's that result in R_nk=0
0192 
0193   // could add an option here to load the munks from file if the dimensions match.
0194 
0195   for (int n = 0; n < NumberOfOrders; n++)  //  !!!  Off by one from Rossegger convention  !!!
0196   {
0197     if (verbosity)
0198     {
0199       std::cout << "Filling Mu[" << n << "][k]..." << std::endl;
0200     }
0201     double x = localepsilon;
0202     for (int k = 0; k < NumberOfOrders; k++)
0203     {
0204       x = FindNextZero(x, localepsilon, n, &Rossegger::Rnk_for_zeroes);
0205       Munk[n][k] = x;
0206       x += localepsilon;
0207       if (verbosity > 0)
0208       {
0209         std::cout << std::format("Mu[{}][{}]={:E}", n, k, Munk[n][k]) << std::endl;
0210     std::cout << std::format("adjacent values are Rnk[mu-localepsilon]={:E}\tRnk[mu+localepsilon]={:E}",
0211                  Rnk_for_zeroes(n, x - localepsilon),
0212                  Rnk_for_zeroes(n, x + localepsilon))
0213           << std::endl;
0214         if (verbosity > 100)
0215         {
0216           std::cout << "values of argument to limu and kimu are " << ((n + 1) * pi / L * a)
0217                     << " and " << ((n + 1) * pi / L * b) << std::endl;
0218         }
0219       }
0220     }
0221   }
0222 
0223   //  Now fill in the N2nk array...
0224   for (int n = 0; n < NumberOfOrders; n++)
0225   {
0226     for (int k = 0; k < NumberOfOrders; k++)
0227     {
0228       //  Rossegger Equation 5.48
0229       //  Integral of R_nk(r)*R_ns(r) dr/r= delta_ks*N2nk
0230       //  note that unlike N2mn, there is no convenient shortcut here.
0231       double integral = 0.0;
0232       double step = 0.001;
0233 
0234       for (double r = a; r < b; r += step)
0235       {
0236         double rnkval = Rnk_(n, k, r);  // must used un-optimized.  We don't have those values yet...
0237 
0238         integral += rnkval * rnkval / r * step;
0239       }
0240       if (verbosity > 1)
0241       {
0242         std::cout << " Int: " << integral << std::endl;
0243       }
0244       N2nk[n][k] = integral;
0245       if (verbosity > 1)
0246       {
0247         std::cout << "n: " << n << " k: " << k << " N2nk[n][k]: " << N2nk[n][k];
0248       }
0249     }
0250   }
0251 
0252   std::cout << "Done." << std::endl;
0253   return;
0254 }
0255 
0256 bool Rossegger::CheckZeroes(double localepsilon)
0257 {
0258   // confirm that the tabulated zeroes are all zeroes of their respective functions:
0259   double result;
0260   for (int m = 0; m < NumberOfOrders; m++)
0261   {
0262     for (int n = 0; n < NumberOfOrders; n++)
0263     {  //  !!!  Off by one from Rossegger convention  !!!
0264       result = Rmn_for_zeroes(m, Betamn[m][n] * b);
0265       if (abs(result) > localepsilon)
0266       {
0267         std::cout << std::format("(m={:d},n={:d}) Jm(x)Ym(lx)-Jm(lx)Ym(x) = {:.6f} for x=b*{:.6f}",
0268                  m, n, result, Betamn[m][n])
0269           << std::endl;
0270         return false;
0271       }
0272     }
0273   }
0274 
0275   // R_nk(mu_nk;a,b)=Limu(Beta_n*a)Kimu(Beta_n*b)-Kimu(Beta_n*a)Limu(Beta_n*b)=0
0276   for (int n = 0; n < NumberOfOrders; n++)
0277   {
0278     for (int k = 0; k < NumberOfOrders; k++)
0279     {  //  !!!  Off by one from Rossegger convention  !!!
0280       result = Rnk_for_zeroes(n, Munk[n][k]);
0281       if (abs(result) > localepsilon * 100)
0282       {
0283         std::cout << std::format("(n={:d},k={:d}) limu(npi*a/L)kimu(npi*b/L)-kimu(npi*a/L)kimu(npi*b/L) = {:.6f} (>eps*100) for mu={:.6f}",
0284                  n, k, result, Munk[n][k])
0285           << std::endl;
0286         return false;
0287       }
0288     }
0289   }
0290 
0291   return true;
0292 }
0293 
0294 void Rossegger::PrecalcFreeConstants()
0295 {  // Routine used to fill the arrays of other values used repeatedly
0296   // these constants depend only on the geometry of the detector
0297   std::cout << "Precalcing " << (3 * NumberOfOrders + 5 * NumberOfOrders * NumberOfOrders)
0298             << "geometric constants" << std::endl;
0299   for (int n = 0; n < NumberOfOrders; n++)
0300   {
0301     BetaN[n] = (n + 1) * pi / L;  //  BetaN=(n+1)*pi/L as used in eg 5.32, .46
0302     BetaN_a[n] = BetaN[n] * a;    //  BetaN*a as used in eg 5.32, .46
0303     BetaN_b[n] = BetaN[n] * b;    //  BetaN*b as used in eg 5.32, .46
0304     for (int m = 0; m < NumberOfOrders; m++)
0305     {
0306       km_BetaN_a[m][n] = kn(m, BetaN_a[n]);                                                                                                                      // kn(m,BetaN*a) as used in Rossegger 5.32
0307       im_BetaN_a[m][n] = in(m, BetaN_a[n]);                                                                                                                      // in(m,BetaN*a) as used in Rossegger 5.32
0308       km_BetaN_b[m][n] = kn(m, BetaN_b[n]);                                                                                                                      // kn(m,BetaN*b) as used in Rossegger 5.33
0309       im_BetaN_b[m][n] = in(m, BetaN_b[n]);                                                                                                                      // in(m,BetaN*b) as used in Rossegger 5.33
0310       bessel_denominator[m][n] = TMath::BesselI(m, BetaN_a[n]) * TMath::BesselK(m, BetaN_b[n]) - TMath::BesselI(m, BetaN_b[n]) * TMath::BesselK(m, BetaN_a[n]);  // TMath::BesselI(m,BetaN*a)*TMath::BesselK(m,BetaN*b)-TMath::BesselI(m,BetaN*b)*TMath::BesselK(m,BetaN*a) as in Rossegger 5.65
0311     }
0312   }
0313   return;
0314 }
0315 void Rossegger::PrecalcDerivedConstants()
0316 {  // Routine used to fill the arrays of other values used repeatedly
0317    // these constants depend on geometry and the zeroes of special functions
0318   std::cout << "Precalcing " << (6 * NumberOfOrders * NumberOfOrders)
0319             << " geometric constants" << std::endl;
0320 
0321   for (int n = 0; n < NumberOfOrders; n++)
0322   {
0323     for (int m = 0; m < NumberOfOrders; m++)
0324     {
0325       ym_Betamn_a[m][n] = yn(m, Betamn[m][n] * a);   // yn(m,Betamn[m][n]*a) as used in Rossegger 5.11
0326       jm_Betamn_a[m][n] = jn(m, Betamn[m][n] * a);   // jn(m,Betamn[m][n]*a) as used in Rossegger 5.11
0327       sinh_Betamn_L[m][n] = sinh(Betamn[m][n] * L);  // sinh(Betamn[m][n]*L)  as in Rossegger 5.64
0328     }
0329     for (int k = 0; k < NumberOfOrders; k++)
0330     {
0331       liMunk_BetaN_a[n][k] = limu(Munk[n][k], BetaN_a[n]);  // limu(Munk[n][k],BetaN*a) as used in Rossegger 5.45
0332       kiMunk_BetaN_a[n][k] = kimu(Munk[n][k], BetaN_a[n]);  // kimu(Munk[n][k],BetaN*a) as used in Rossegger 5.45
0333       sinh_pi_Munk[n][k] = sinh(pi * Munk[n][k]);           // sinh(pi*Munk[n][k]) as in Rossegger 5.66
0334     }
0335   }
0336   return;
0337 }
0338 
0339 double Rossegger::Limu(double mu, double x)
0340 {
0341   // defined in Rossegger eqn 5.44, also a canonical 'satisfactory companion' to Kimu.
0342   // could use Griddit?
0343   //   Rossegger Equation 5.45
0344   //        Rnk(r) = Limu_nk (BetaN a) Kimu_nk (BetaN r) - Kimu_nk(BetaN a) Limu_nk (BetaN r)
0345 
0346   int IFAC = 1;
0347   double A = mu;
0348   double DLI = 0;
0349   double DERR = 0;
0350   int IERRO = 0;
0351 
0352   double X = x;
0353   dlia_(&IFAC, &X, &A, &DLI, &DERR, &IERRO);
0354   return DLI;
0355 }
0356 double Rossegger::Kimu(double mu, double x)
0357 {
0358   int IFAC = 1;
0359   double A = mu;
0360   double DKI = 0;
0361   double DERR = 0;
0362   int IERRO = 0;
0363 
0364   double X = x;
0365   dkia_(&IFAC, &X, &A, &DKI, &DERR, &IERRO);
0366   return DKI;
0367 }
0368 
0369 double Rossegger::Rmn_for_zeroes(int m, double x) // NOLINT(readability-make-member-function-const)
0370 {
0371   double lx = a * x / b;
0372   //  Rossegger Equation 5.12:
0373 
0374   return jn(m, x) * yn(m, lx) - jn(m, lx) * yn(m, x);
0375 }
0376 
0377 double Rossegger::Rmn(int m, int n, double r)
0378 {
0379   if (verbosity > 100)
0380   {
0381     std::cout << "Determine Rmn(" << m << "," << n << "," << r << ") = ";
0382   }
0383 
0384   //  Check input arguments for sanity...
0385   int error = 0;
0386   if (m < 0 || m >= NumberOfOrders)
0387   {
0388     error = 1;
0389   }
0390   if (n < 0 || n >= NumberOfOrders)
0391   {
0392     error = 1;
0393   }
0394   if (r < a || r > b)
0395   {
0396     error = 1;
0397   }
0398   if (error)
0399   {
0400     std::cout << "Invalid arguments Rmn(" << m << "," << n << "," << r << ")" << std::endl;
0401     ;
0402     return 0;
0403   }
0404 
0405   //  Calculate the function using C-libraries from boost
0406   //  Rossegger Equation 5.11:
0407   //         Rmn(r) = Ym(Beta_mn a)*Jm(Beta_mn r) - Jm(Beta_mn a)*Ym(Beta_mn r)
0408   double R = 0;
0409   R = ym_Betamn_a[m][n] * jn(m, Betamn[m][n] * r) - jm_Betamn_a[m][n] * yn(m, Betamn[m][n] * r);
0410 
0411   if (verbosity > 100)
0412   {
0413     std::cout << R << std::endl;
0414   }
0415   return R;
0416 }
0417 
0418 double Rossegger::Rmn_(int m, int n, double r)
0419 {
0420   if (verbosity > 100)
0421   {
0422     std::cout << "Determine Rmn(" << m << "," << n << "," << r << ") = ";
0423   }
0424 
0425   //  Check input arguments for sanity...
0426   int error = 0;
0427   if (m < 0 || m >= NumberOfOrders)
0428   {
0429     error = 1;
0430   }
0431   if (n < 0 || n >= NumberOfOrders)
0432   {
0433     error = 1;
0434   }
0435   if (r < a || r > b)
0436   {
0437     error = 1;
0438   }
0439   if (error)
0440   {
0441     std::cout << "Invalid arguments Rmn(" << m << "," << n << "," << r << ")" << std::endl;
0442     ;
0443     return 0;
0444   }
0445 
0446   //  Calculate the function using C-libraries from boost
0447   //  Rossegger Equation 5.11:
0448   //         Rmn(r) = Ym(Beta_mn a)*Jm(Beta_mn r) - Jm(Beta_mn a)*Ym(Beta_mn r)
0449   double R = 0;
0450   R = yn(m, Betamn[m][n] * a) * jn(m, Betamn[m][n] * r) - jn(m, Betamn[m][n] * a) * yn(m, Betamn[m][n] * r);
0451 
0452   if (verbosity > 100)
0453   {
0454     std::cout << R << std::endl;
0455   }
0456   return R;
0457 }
0458 
0459 double Rossegger::Rmn1(int m, int n, double r)
0460 {
0461   //  Check input arguments for sanity...
0462   int error = 0;
0463   if (m < 0 || m >= NumberOfOrders)
0464   {
0465     error = 1;
0466   }
0467   if (n < 0 || n >= NumberOfOrders)
0468   {
0469     error = 1;
0470   }
0471   if (r < a || r > b)
0472   {
0473     error = 1;
0474   }
0475   if (error)
0476   {
0477     std::cout << "Invalid arguments Rmn1(" << m << "," << n << "," << r << ")" << std::endl;
0478     ;
0479     return 0;
0480   }
0481 
0482   //  Calculate using the TMath functions from root.
0483   //  Rossegger Equation 5.32
0484   //         Rmn1(r) = Km(BetaN a)Im(BetaN r) - Im(BetaN a) Km(BetaN r)
0485   double R = 0;
0486   R = km_BetaN_a[m][n] * in(m, BetaN[n] * r) - im_BetaN_a[m][n] * kn(m, BetaN[n] * r);
0487 
0488   return R;
0489 }
0490 
0491 double Rossegger::Rmn1_(int m, int n, double r) // NOLINT(readability-make-member-function-const)
0492 {
0493   //  Check input arguments for sanity...
0494   int error = 0;
0495   if (m < 0 || m >= NumberOfOrders)
0496   {
0497     error = 1;
0498   }
0499   if (n < 0 || n >= NumberOfOrders)
0500   {
0501     error = 1;
0502   }
0503   if (r < a || r > b)
0504   {
0505     error = 1;
0506   }
0507   if (error)
0508   {
0509     std::cout << "Invalid arguments Rmn1(" << m << "," << n << "," << r << ")" << std::endl;
0510     ;
0511     return 0;
0512   }
0513 
0514   //  Calculate using the TMath functions from root.
0515   //  Rossegger Equation 5.32
0516   //         Rmn1(r) = Km(BetaN a)Im(BetaN r) - Im(BetaN a) Km(BetaN r)
0517   double R = 0;
0518   double BetaN_ = (n + 1) * pi / L;
0519   R = kn(m, BetaN_ * a) * in(m, BetaN_ * r) - in(m, BetaN_ * a) * kn(m, BetaN_ * r);
0520 
0521   return R;
0522 }
0523 
0524 double Rossegger::Rmn2(int m, int n, double r)
0525 {
0526   //  Check input arguments for sanity...
0527   int error = 0;
0528   if (m < 0 || m >= NumberOfOrders)
0529   {
0530     error = 1;
0531   }
0532   if (n < 0 || n >= NumberOfOrders)
0533   {
0534     error = 1;
0535   }
0536   if (r < a || r > b)
0537   {
0538     error = 1;
0539   }
0540   if (error)
0541   {
0542     std::cout << "Invalid arguments Rmn2(" << m << "," << n << "," << r << ")" << std::endl;
0543     ;
0544     return 0;
0545   }
0546 
0547   //  Calculate using the TMath functions from root.
0548   //  Rossegger Equation 5.33
0549   //         Rmn2(r) = Km(BetaN b)Im(BetaN r) - Im(BetaN b) Km(BetaN r)
0550   double R = 0;
0551   R = km_BetaN_b[m][n] * in(m, BetaN[n] * r) - im_BetaN_b[m][n] * kn(m, BetaN[n] * r);
0552 
0553   return R;
0554 }
0555 
0556 double Rossegger::Rmn2_(int m, int n, double r) // NOLINT(readability-make-member-function-const)
0557 {
0558   //  Check input arguments for sanity...
0559   int error = 0;
0560   if (m < 0 || m >= NumberOfOrders)
0561   {
0562     error = 1;
0563   }
0564   if (n < 0 || n >= NumberOfOrders)
0565   {
0566     error = 1;
0567   }
0568   if (r < a || r > b)
0569   {
0570     error = 1;
0571   }
0572   if (error)
0573   {
0574     std::cout << "Invalid arguments Rmn2(" << m << "," << n << "," << r << ")" << std::endl;
0575     ;
0576     return 0;
0577   }
0578 
0579   //  Calculate using the TMath functions from root.
0580   //  Rossegger Equation 5.33
0581   //         Rmn2(r) = Km(BetaN b)Im(BetaN r) - Im(BetaN b) Km(BetaN r)
0582   double R = 0;
0583   double BetaN_ = (n + 1) * pi / L;
0584   R = kn(m, BetaN_ * b) * in(m, BetaN_ * r) - in(m, BetaN_ * b) * kn(m, BetaN_ * r);
0585 
0586   return R;
0587 }
0588 
0589 double Rossegger::RPrime(int m, int n, double ref, double r)
0590 {
0591   //  Check input arguments for sanity...
0592   int error = 0;
0593   if (m < 0 || m >= NumberOfOrders)
0594   {
0595     error = 1;
0596   }
0597   if (n < 0 || n >= NumberOfOrders)
0598   {
0599     error = 1;
0600   }
0601   if (ref < a || ref > b)
0602   {
0603     error = 1;
0604   }
0605   if (r < a || r > b)
0606   {
0607     error = 1;
0608   }
0609   if (error)
0610   {
0611     std::cout << "Invalid arguments RPrime(" << m << "," << n << "," << ref << "," << r << ")" << std::endl;
0612     ;
0613     return 0;
0614   }
0615 
0616   double R = 0;
0617   //  Calculate using the TMath functions from root.
0618   //  Rossegger Equation 5.65
0619   //         Rmn2(ref,r) = BetaN/2* [   Km(BetaN ref) {Im-1(BetaN r) + Im+1(BetaN r)}
0620   //                                  - Im(BetaN ref) {Km-1(BetaN r) + Km+1(BetaN r)}  ]
0621   //  NOTE:  K-m(z) = Km(z) and I-m(z) = Im(z)... though boost handles negative orders.
0622   //
0623   // with: s -> ref,  t -> r,
0624   double BetaN_ = BetaN[n];
0625   double term1 = kn(m, BetaN_ * ref) * (in(m - 1, BetaN_ * r) + in(m + 1, BetaN_ * r));
0626   double term2 = in(m, BetaN_ * ref) * (kn(m - 1, BetaN_ * r) + kn(m + 1, BetaN_ * r));
0627   R = BetaN_ / 2.0 * (term1 + term2);
0628 
0629   return R;
0630 }
0631 
0632 double Rossegger::RPrime_(int m, int n, double ref, double r) // NOLINT(readability-make-member-function-const)
0633 {
0634   //  Check input arguments for sanity...
0635   int error = 0;
0636   if (m < 0 || m >= NumberOfOrders)
0637   {
0638     error = 1;
0639   }
0640   if (n < 0 || n >= NumberOfOrders)
0641   {
0642     error = 1;
0643   }
0644   if (ref < a || ref > b)
0645   {
0646     error = 1;
0647   }
0648   if (r < a || r > b)
0649   {
0650     error = 1;
0651   }
0652   if (error)
0653   {
0654     std::cout << "Invalid arguments RPrime(" << m << "," << n << "," << ref << "," << r << ")" << std::endl;
0655     ;
0656     return 0;
0657   }
0658 
0659   double R = 0;
0660   //  Rossegger Equation 5.65
0661   //         Rmn2(ref,r) = BetaN/2* [   Km(BetaN ref) {Im-1(BetaN r) + Im+1(BetaN r)}
0662   //                                  - Im(BetaN ref) {Km-1(BetaN r) + Km+1(BetaN r)}  ]
0663   //  NOTE:  K-m(z) = Km(z) and I-m(z) = Im(z)... though boost handles negative orders.
0664   //
0665   // with: s -> ref,  t -> r,
0666   double BetaN_ = (n + 1) * pi / L;
0667   double term1 = kn(m, BetaN_ * ref) * (in(m - 1, BetaN_ * r) + in(m + 1, BetaN_ * r));
0668   double term2 = in(m, BetaN_ * ref) * (kn(m - 1, BetaN_ * r) + kn(m + 1, BetaN_ * r));
0669   R = BetaN_ / 2.0 * (term1 + term2);
0670 
0671   return R;
0672 }
0673 
0674 double Rossegger::Rnk_for_zeroes(int n, double mu)
0675 {
0676   // unlike Rossegger, we count 'k' and 'n' from zero.
0677   if (verbosity > 10)
0678   {
0679     std::cout << "Rnk_for_zeroes called with n=" << n << ",mu=" << mu << std::endl;
0680   }
0681   double betana = BetaN_a[n];
0682   double betanb = BetaN_b[n];
0683   //  Rossegger Equation 5.46
0684   //       Rnk(r) = Limu_nk (BetaN a) Kimu_nk (BetaN b) - Kimu_nk(BetaN a) Limu_nk (BetaN b)
0685 
0686   return limu(mu, betana) * kimu(mu, betanb) - kimu(mu, betana) * limu(mu, betanb);
0687 }
0688 
0689 double Rossegger::Rnk_for_zeroes_(int n, double mu)
0690 {
0691   // unlike Rossegger, we count 'k' and 'n' from zero.
0692   if (verbosity > 10)
0693   {
0694     std::cout << "Rnk_for_zeroes called with n=" << n << ",mu=" << mu << std::endl;
0695   }
0696   double BetaN_ = (n + 1) * pi / L;  // this is defined in the paragraph before 5.46
0697   double betana = BetaN_ * a;
0698   double betanb = BetaN_ * b;
0699   //  Rossegger Equation 5.46
0700   //       Rnk(r) = Limu_nk (BetaN a) Kimu_nk (BetaN b) - Kimu_nk(BetaN a) Limu_nk (BetaN b)
0701 
0702   return limu(mu, betana) * kimu(mu, betanb) - kimu(mu, betana) * limu(mu, betanb);
0703 }
0704 
0705 double Rossegger::Rnk(int n, int k, double r)
0706 {
0707   //  Check input arguments for sanity...
0708   int error = 0;
0709   if (n < 0 || n >= NumberOfOrders)
0710   {
0711     error = 1;
0712   }
0713   if (k < 0 || k >= NumberOfOrders)
0714   {
0715     error = 1;
0716   }
0717   if (r < a || r > b)
0718   {
0719     error = 1;
0720   }
0721   if (error)
0722   {
0723     std::cout << "Invalid arguments Rnk(" << n << "," << k << "," << r << ")" << std::endl;
0724     ;
0725     return 0;
0726   }
0727   //  Rossegger Equation 5.45
0728   //       Rnk(r) = Limu_nk (BetaN a) Kimu_nk (BetaN r) - Kimu_nk(BetaN a) Limu_nk (BetaN r)
0729 
0730   return liMunk_BetaN_a[n][k] * kimu(Munk[n][k], BetaN[n] * r) - kiMunk_BetaN_a[n][k] * limu(Munk[n][k], BetaN[n] * r);
0731 }
0732 
0733 double Rossegger::Rnk_(int n, int k, double r)
0734 {
0735   //  Check input arguments for sanity...
0736   int error = 0;
0737   if (n < 0 || n >= NumberOfOrders)
0738   {
0739     error = 1;
0740   }
0741   if (k < 0 || k >= NumberOfOrders)
0742   {
0743     error = 1;
0744   }
0745   if (r < a || r > b)
0746   {
0747     error = 1;
0748   }
0749   if (error)
0750   {
0751     std::cout << "Invalid arguments Rnk(" << n << "," << k << "," << r << ")" << std::endl;
0752     ;
0753     return 0;
0754   }
0755   double BetaN_ = (n + 1) * pi / L;
0756   //  Rossegger Equation 5.45
0757   //       Rnk(r) = Limu_nk (BetaN a) Kimu_nk (BetaN r) - Kimu_nk(BetaN a) Limu_nk (BetaN r)
0758 
0759   return limu(Munk[n][k], BetaN_ * a) * kimu(Munk[n][k], BetaN_ * r) - kimu(Munk[n][k], BetaN_ * a) * limu(Munk[n][k], BetaN_ * r);
0760 }
0761 
0762 double Rossegger::Ez(double r, double phi, double z, double r1, double phi1, double z1)
0763 {
0764   // rcc streamlined Ez
0765   int error = 0;
0766   if (r < a || r > b)
0767   {
0768     error = 1;
0769   }
0770   if (phi < 0 || phi > 2 * pi)
0771   {
0772     error = 1;
0773   }
0774   if (z < 0 || z > L)
0775   {
0776     error = 1;
0777   }
0778   if (r1 < a || r1 > b)
0779   {
0780     error = 1;
0781   }
0782   if (phi1 < 0 || phi1 > 2 * pi)
0783   {
0784     error = 1;
0785   }
0786   if (z1 < 0 || z1 > L)
0787   {
0788     error = 1;
0789   }
0790   if (error)
0791   {
0792     std::cout << "Invalid arguments Ez(";
0793     std::cout << r << ",";
0794     std::cout << phi << ",";
0795     std::cout << z << ",";
0796     std::cout << r1 << ",";
0797     std::cout << phi1 << ",";
0798     std::cout << z1;
0799     std::cout << ")" << std::endl;
0800     ;
0801     return 0;
0802   }
0803   // Rossegger Equation 5.64
0804   double G = 0;
0805   for (int m = 0; m < NumberOfOrders; m++)
0806   {
0807     if (verbosity > 10)
0808     {
0809       std::cout << std::endl
0810                 << m;
0811     }
0812     for (int n = 0; n < NumberOfOrders; n++)
0813     {
0814       if (verbosity > 10)
0815       {
0816         std::cout << " " << n;
0817       }
0818       double term = 1;  // unitless
0819       if (verbosity > 10)
0820       {
0821         std::cout << " " << term;
0822       }
0823       term *= (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));  // unitless
0824       if (verbosity > 10)
0825       {
0826         std::cout << " " << term;
0827       }
0828       term *= Rmn(m, n, r) * Rmn(m, n, r1) / N2mn[m][n];  // units of 1/[L]^2
0829       if (verbosity > 10)
0830       {
0831         std::cout << " " << term;
0832       }
0833       if (z < z1)
0834       {
0835         term *= cosh(Betamn[m][n] * z) * sinh(Betamn[m][n] * (L - z1)) / sinh_Betamn_L[m][n];  // unitless
0836       }
0837       else
0838       {
0839         term *= -cosh(Betamn[m][n] * (L - z)) * sinh(Betamn[m][n] * z1) / sinh_Betamn_L[m][n];  // unitless
0840       }
0841       if (verbosity > 10)
0842       {
0843         std::cout << " " << term;
0844       }
0845       G += term;
0846       if (verbosity > 10)
0847       {
0848         std::cout << " " << term << " " << G << std::endl;
0849       }
0850     }
0851   }
0852 
0853   G = G / (2.0 * pi);
0854   if (verbosity)
0855   {
0856     std::cout << "Ez = " << G << std::endl;
0857   }
0858 
0859   return G;
0860 }
0861 
0862 double Rossegger::Ez_(double r, double phi, double z, double r1, double phi1, double z1)
0863 {
0864   // if(fByFile && fabs(r-r1)>MinimumDR && fabs(z-z1)>MinimumDZ) return ByFileEZ(r,phi,z,r1,phi1,z1);
0865   //   Check input arguments for sanity...
0866   int error = 0;
0867   if (r < a || r > b)
0868   {
0869     error = 1;
0870   }
0871   if (phi < 0 || phi > 2 * pi)
0872   {
0873     error = 1;
0874   }
0875   if (z < 0 || z > L)
0876   {
0877     error = 1;
0878   }
0879   if (r1 < a || r1 > b)
0880   {
0881     error = 1;
0882   }
0883   if (phi1 < 0 || phi1 > 2 * pi)
0884   {
0885     error = 1;
0886   }
0887   if (z1 < 0 || z1 > L)
0888   {
0889     error = 1;
0890   }
0891   if (error)
0892   {
0893     std::cout << "Invalid arguments Ez(";
0894     std::cout << r << ",";
0895     std::cout << phi << ",";
0896     std::cout << z << ",";
0897     std::cout << r1 << ",";
0898     std::cout << phi1 << ",";
0899     std::cout << z1;
0900     std::cout << ")" << std::endl;
0901     ;
0902     return 0;
0903   }
0904 
0905   double G = 0;
0906   for (int m = 0; m < NumberOfOrders; m++)
0907   {
0908     if (verbosity > 10)
0909     {
0910       std::cout << std::endl
0911                 << m;
0912     }
0913     for (int n = 0; n < NumberOfOrders; n++)
0914     {
0915       if (verbosity > 10)
0916       {
0917         std::cout << " " << n;
0918       }
0919       double term = 1 / (2.0 * pi);
0920       if (verbosity > 10)
0921       {
0922         std::cout << " " << term;
0923       }
0924       term *= (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));  // unitless
0925       if (verbosity > 10)
0926       {
0927         std::cout << " " << term;
0928       }
0929       term *= Rmn(m, n, r) * Rmn(m, n, r1) / N2mn[m][n];  // units of 1/[L]^2
0930       if (verbosity > 10)
0931       {
0932         std::cout << " " << term;
0933       }
0934       if (z < z1)
0935       {
0936         term *= cosh(Betamn[m][n] * z) * sinh(Betamn[m][n] * (L - z1)) / sinh(Betamn[m][n] * L);
0937       }
0938       else
0939       {
0940         term *= -cosh(Betamn[m][n] * (L - z)) * sinh(Betamn[m][n] * z1) / sinh(Betamn[m][n] * L);
0941         ;
0942       }
0943       if (verbosity > 10)
0944       {
0945         std::cout << " " << term;
0946       }
0947       G += term;
0948       if (verbosity > 10)
0949       {
0950         std::cout << " " << term << " " << G << std::endl;
0951       }
0952     }
0953   }
0954   if (verbosity)
0955   {
0956     std::cout << "Ez = " << G << std::endl;
0957   }
0958 
0959   return G;
0960 }
0961 
0962 double Rossegger::Er(double r, double phi, double z, double r1, double phi1, double z1)
0963 {
0964   // rcc streamlined Er
0965 
0966   // as in Rossegger 5.65
0967   // field at r, phi, z due to unit charge at r1, phi1, z1;
0968   // if(fByFile && fabs(r-r1)>MinimumDR && fabs(z-z1)>MinimumDZ) return ByFileER(r,phi,z,r1,phi1,z1);
0969   //   Check input arguments for sanity...
0970   int error = 0;
0971   if (r < a || r > b)
0972   {
0973     error = 1;
0974   }
0975   if (phi < 0 || phi > 2 * pi)
0976   {
0977     error = 1;
0978   }
0979   if (z < 0 || z > L)
0980   {
0981     error = 1;
0982   }
0983   if (r1 < a || r1 > b)
0984   {
0985     error = 1;
0986   }
0987   if (phi1 < 0 || phi1 > 2 * pi)
0988   {
0989     error = 1;
0990   }
0991   if (z1 < 0 || z1 > L)
0992   {
0993     error = 1;
0994   }
0995   if (error)
0996   {
0997     std::cout << "Invalid arguments Er(";
0998     std::cout << r << ",";
0999     std::cout << phi << ",";
1000     std::cout << z << ",";
1001     std::cout << r1 << ",";
1002     std::cout << phi1 << ",";
1003     std::cout << z1;
1004     std::cout << ")" << std::endl;
1005     ;
1006     return 0;
1007   }
1008 
1009   double part = 0;
1010   double G = 0;
1011   for (int m = 0; m < NumberOfOrders; m++)
1012   {
1013     for (int n = 0; n < NumberOfOrders; n++)
1014     {
1015       double term = 1;
1016       part = (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));  // unitless
1017       if (verbosity > 10)
1018       {
1019         std::cout << "(2 - ((m==0)?1:0))*cos(m*(phi-phi1)); = " << part << std::endl;
1020       }
1021       term *= part;
1022       part = sin(BetaN[n] * z) * sin(BetaN[n] * z1);  // unitless
1023       if (verbosity > 10)
1024       {
1025         std::cout << "sin(BetaN[n]*z)*sin(BetaN[n]*z1); = " << part << std::endl;
1026       }
1027       term *= part;
1028 
1029       if (r < r1)
1030       {
1031         term *= RPrime(m, n, a, r) * Rmn2(m, n, r1);  // units of 1/[L]
1032       }
1033       else
1034       {
1035         term *= Rmn1(m, n, r1) * RPrime(m, n, b, r);  // units of 1/[L]
1036       }
1037       term /= bessel_denominator[m][n];  // unitless
1038       G += term;
1039     }
1040   }
1041 
1042   G = G / (L * pi);  // units of 1/[L] -- net is 1/[L]^2
1043   if (verbosity)
1044   {
1045     std::cout << "Er = " << G << std::endl;
1046   }
1047 
1048   return G;
1049 }
1050 
1051 double Rossegger::Er_(double r, double phi, double z, double r1, double phi1, double z1)
1052 {
1053   // doesn't take advantage of precalcs.
1054   // as in Rossegger 5.65
1055   // field at r, phi, z due to unit charge at r1, phi1, z1;
1056   // if(fByFile && fabs(r-r1)>MinimumDR && fabs(z-z1)>MinimumDZ) return ByFileER(r,phi,z,r1,phi1,z1);
1057   //   Check input arguments for sanity...
1058   int error = 0;
1059   if (r < a || r > b)
1060   {
1061     error = 1;
1062   }
1063   if (phi < 0 || phi > 2 * pi)
1064   {
1065     error = 1;
1066   }
1067   if (z < 0 || z > L)
1068   {
1069     error = 1;
1070   }
1071   if (r1 < a || r1 > b)
1072   {
1073     error = 1;
1074   }
1075   if (phi1 < 0 || phi1 > 2 * pi)
1076   {
1077     error = 1;
1078   }
1079   if (z1 < 0 || z1 > L)
1080   {
1081     error = 1;
1082   }
1083   if (error)
1084   {
1085     std::cout << "Invalid arguments Er(";
1086     std::cout << r << ",";
1087     std::cout << phi << ",";
1088     std::cout << z << ",";
1089     std::cout << r1 << ",";
1090     std::cout << phi1 << ",";
1091     std::cout << z1;
1092     std::cout << ")" << std::endl;
1093     ;
1094     return 0;
1095   }
1096 
1097   double G = 0;
1098   for (int m = 0; m < NumberOfOrders; m++)
1099   {
1100     for (int n = 0; n < NumberOfOrders; n++)
1101     {
1102       double term = 1 / (L * pi);
1103       term *= (2 - ((m == 0) ? 1 : 0)) * cos(m * (phi - phi1));
1104       double BetaN_ = (n + 1) * pi / L;
1105       term *= sin(BetaN_ * z) * sin(BetaN_ * z1);
1106       if (r < r1)
1107       {
1108         term *= RPrime_(m, n, a, r) * Rmn2_(m, n, r1);
1109       }
1110       else
1111       {
1112         term *= Rmn1_(m, n, r1) * RPrime_(m, n, b, r);
1113       }
1114 
1115       term /= TMath::BesselI(m, BetaN_ * a) * TMath::BesselK(m, BetaN_ * b) - TMath::BesselI(m, BetaN_ * b) * TMath::BesselK(m, BetaN_ * a);
1116 
1117       G += term;
1118     }
1119   }
1120 
1121   if (verbosity)
1122   {
1123     std::cout << "Er = " << G << std::endl;
1124   }
1125 
1126   return G;
1127 }
1128 
1129 double Rossegger::Ephi(double r, double phi, double z, double r1, double phi1, double z1)
1130 {
1131   // rcc streamlined Ephi term
1132   // compute field at rphiz from charge at r1phi1z1
1133   //   Check input arguments for sanity...
1134   int error = 0;
1135   if (r < a || r > b)
1136   {
1137     error = 1;
1138   }
1139   if (phi < 0 || phi > 2 * pi)
1140   {
1141     error = 1;
1142   }
1143   if (z < 0 || z > L)
1144   {
1145     error = 1;
1146   }
1147   if (r1 < a || r1 > b)
1148   {
1149     error = 1;
1150   }
1151   if (phi1 < 0 || phi1 > 2 * pi)
1152   {
1153     error = 1;
1154   }
1155   if (z1 < 0 || z1 > L)
1156   {
1157     error = 1;
1158   }
1159   if (error)
1160   {
1161     std::cout << "Invalid arguments Ephi(";
1162     std::cout << r << ",";
1163     std::cout << phi << ",";
1164     std::cout << z << ",";
1165     std::cout << r1 << ",";
1166     std::cout << phi1 << ",";
1167     std::cout << z1;
1168     std::cout << ")" << std::endl;
1169     ;
1170     return 0;
1171   }
1172 
1173   double G = 0;
1174   // Rossegger Eqn. 5.66:
1175   for (int k = 0; k < NumberOfOrders; k++)
1176   {
1177     for (int n = 0; n < NumberOfOrders; n++)
1178     {
1179       double term = 1;
1180       term *= sin(BetaN[n] * z) * sin(BetaN[n] * z1);     // unitless
1181       term *= Rnk(n, k, r) * Rnk(n, k, r1) / N2nk[n][k];  // unitless?
1182 
1183       // the derivative of cosh(munk(pi-|phi-phi1|)
1184       if (phi > phi1)
1185       {
1186         term *= -sinh(Munk[n][k] * (pi - (phi - phi1)));  // unitless
1187       }
1188       else
1189       {
1190         term *= sinh(Munk[n][k] * (pi - (phi1 - phi)));  // unitless
1191       }
1192       term *= 1 / sinh_pi_Munk[n][k];  // unitless
1193       G += term;
1194     }
1195   }
1196 
1197   G = G / (L * r);  // units of 1/[L]^2.  r comes from the phi term in cylindrical gradient expression.
1198   if (verbosity)
1199   {
1200     std::cout << "Ephi = " << G << std::endl;
1201   }
1202 
1203   return G;
1204 }
1205 
1206 double Rossegger::Ephi_(double r, double phi, double z, double r1, double phi1, double z1)
1207 {
1208   // compute field at rphiz from charge at r1phi1z1
1209   //   Check input arguments for sanity...
1210   int error = 0;
1211   if (r < a || r > b)
1212   {
1213     error = 1;
1214   }
1215   if (phi < 0 || phi > 2 * pi)
1216   {
1217     error = 1;
1218   }
1219   if (z < 0 || z > L)
1220   {
1221     error = 1;
1222   }
1223   if (r1 < a || r1 > b)
1224   {
1225     error = 1;
1226   }
1227   if (phi1 < 0 || phi1 > 2 * pi)
1228   {
1229     error = 1;
1230   }
1231   if (z1 < 0 || z1 > L)
1232   {
1233     error = 1;
1234   }
1235   if (error)
1236   {
1237     std::cout << "Invalid arguments Ephi(";
1238     std::cout << r << ",";
1239     std::cout << phi << ",";
1240     std::cout << z << ",";
1241     std::cout << r1 << ",";
1242     std::cout << phi1 << ",";
1243     std::cout << z1;
1244     std::cout << ")" << std::endl;
1245     ;
1246     return 0;
1247   }
1248 
1249   verbosity = 1;
1250   double G = 0;
1251   // Rossegger Eqn. 5.66:
1252   for (int k = 0; k < NumberOfOrders; k++)  // off by one from Rossegger convention!
1253   {
1254     if (verbosity)
1255     {
1256       std::cout << "\nk=" << k;
1257     }
1258     for (int n = 0; n < NumberOfOrders; n++)  // off by one from Rossegger convention!
1259     {
1260       if (verbosity)
1261       {
1262         std::cout << " n=" << n;
1263       }
1264       double BetaN_ = (n + 1) * pi / L;
1265       double term = 1 / (L * r);
1266       if (verbosity)
1267       {
1268         std::cout << " 1/L=" << term;
1269       }
1270       term *= sin(BetaN_ * z) * sin(BetaN_ * z1);
1271       if (verbosity)
1272       {
1273         std::cout << " *sinsin=" << term;
1274       }
1275       term *= Rnk_(n, k, r) * Rnk_(n, k, r1);
1276       if (verbosity)
1277       {
1278         std::cout << " *rnkrnk=" << term;
1279       }
1280       term /= N2nk[n][k];
1281       if (verbosity)
1282       {
1283         std::cout << " */nnknnk=" << term;
1284       }
1285 
1286       // the derivative of cosh(munk(pi-|phi-phi1|)
1287       if (phi > phi1)
1288       {
1289         term *= -sinh(Munk[n][k] * (pi - (phi - phi1)));
1290         // term *=  Munk[n][k]*sinh(Munk[n][k]*pi*(phi1-phi));
1291         // this originally has a factor of Munk in front, but that cancels with one in the denominator
1292       }
1293       else
1294       {
1295         term *= sinh(Munk[n][k] * (pi - (phi1 - phi)));
1296         // term *=  -Munk[n][k]*sinh(Munk[n][k]*pi*(phi-phi1));
1297         // this originally has a factor of Munk in front, but that cancels with one in the denominator
1298       }
1299       if (verbosity)
1300       {
1301         std::cout << " *sinh(mu*pi-phi-phi)=" << term;
1302       }
1303       term *= 1 / (sinh(pi * Munk[n][k]));
1304       // term *= 1/(Munk[n][k]*sinh(pi*Munk[n][k]));
1305       // this originally has a factor of Munk in front, but that cancels with one in the numerator
1306       G += term;
1307       if (verbosity)
1308       {
1309         std::cout << "  /sinh=" << term << " G=" << G << std::endl;
1310       }
1311     }
1312   }
1313   if (verbosity)
1314   {
1315     std::cout << "Ephi = " << G << std::endl;
1316   }
1317   verbosity = 0;
1318 
1319   return G;
1320 }
1321 
1322 void Rossegger::SaveZeroes(const std::string &destfile)
1323 {
1324   TFile *output = TFile::Open(destfile.c_str(), "RECREATE");
1325   output->cd();
1326 
1327   TTree *tInfo = new TTree("info", "Mu[n][k] values");
1328   int ord = NumberOfOrders;
1329   tInfo->Branch("order", &ord);
1330   tInfo->Branch("epsilon", &epsilon);
1331   tInfo->Fill();
1332 
1333   int n;
1334   int k;
1335   int m;
1336   double munk;
1337   double betamn;
1338   double n2nk;
1339   double n2mn;
1340   TTree *tmunk = new TTree("munk", "Mu[n][k] values");
1341   tmunk->Branch("n", &n);
1342   tmunk->Branch("k", &k);
1343   tmunk->Branch("munk", &munk);
1344   tmunk->Branch("n2nk", &n2nk);
1345   for (n = 0; n < ord; n++)
1346   {
1347     for (k = 0; k < ord; k++)
1348     {
1349       munk = Munk[n][k];
1350       n2nk = N2nk[n][k];
1351       tmunk->Fill();
1352     }
1353   }
1354 
1355   TTree *tbetamn = new TTree("betamn", "Beta[m][n] values");
1356   tbetamn->Branch("m", &m);
1357   tbetamn->Branch("n", &n);
1358   tbetamn->Branch("betamn", &betamn);
1359   tbetamn->Branch("n2mn", &n2mn);
1360   for (m = 0; m < ord; m++)
1361   {
1362     for (n = 0; n < ord; n++)
1363     {
1364       betamn = Betamn[m][n];
1365       n2mn = N2mn[m][n];
1366       tbetamn->Fill();
1367     }
1368   }
1369 
1370   tInfo->Write();
1371   tmunk->Write();
1372   tbetamn->Write();
1373   // output->Write();
1374   output->Close();
1375   return;
1376 }
1377 
1378 void Rossegger::LoadZeroes(const std::string &destfile)
1379 {
1380   TFile *f = TFile::Open(destfile.c_str(), "READ");
1381   std::cout << "reading rossegger zeroes from " << destfile << std::endl;
1382   TTree *tInfo = (TTree *) (f->Get("info"));
1383   int ord;
1384   tInfo->SetBranchAddress("order", &ord);
1385   tInfo->SetBranchAddress("epsilon", &epsilon);
1386   tInfo->GetEntry(0);
1387   std::cout << "order=" << ord << ",epsilon=" << epsilon << std::endl;
1388 
1389   int n;
1390   int k;
1391   int m;
1392   double munk;
1393   double betamn;
1394   double n2nk;
1395   double n2mn;
1396   TTree *tmunk = (TTree *) (f->Get("munk"));
1397   tmunk->SetBranchAddress("n", &n);
1398   tmunk->SetBranchAddress("k", &k);
1399   tmunk->SetBranchAddress("munk", &munk);
1400   tmunk->SetBranchAddress("n2nk", &n2nk);
1401   for (int i = 0; i < tmunk->GetEntries(); i++)
1402   {
1403     tmunk->GetEntry(i);
1404     Munk[n][k] = munk;
1405     N2nk[n][k] = n2nk;
1406   }
1407 
1408   TTree *tbetamn = (TTree *) (f->Get("betamn"));
1409   tbetamn->SetBranchAddress("m", &m);
1410   tbetamn->SetBranchAddress("n", &n);
1411   tbetamn->SetBranchAddress("betamn", &betamn);
1412   tbetamn->SetBranchAddress("n2mn", &n2mn);
1413   for (int i = 0; i < tbetamn->GetEntries(); i++)
1414   {
1415     tbetamn->GetEntry(i);
1416     Betamn[m][n] = betamn;
1417     N2mn[m][n] = n2mn;
1418   }
1419 
1420   f->Close();
1421   return;
1422 }