Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:15:40

0001 #ifndef TPCGENERATOR_ROSSEGGER_H
0002 #define TPCGENERATOR_ROSSEGGER_H
0003 
0004 //
0005 //  Hello Space Charge Fans:  (although you should hate space charge)
0006 //
0007 //    This is a code that implements the calculations of space charge contributions
0008 //  In a cylindrical TPC.  It uses the solutions discovered by Stefan Rossegger for
0009 //  the ALICE experiment.  These solutions use various Bessel functions as a series
0010 //  solution to the "point charge in a conducting cylinder" problem imposed by
0011 //  The typical configuration of a TPC found at a collider.
0012 //
0013 //    These calculations have only a single dimension [length] and we shall choose
0014 //  the units of cm throughout the calculations.
0015 //
0016 //                                                TKH
0017 //                                                12-3-2015
0018 //
0019 
0020 #include <cmath>
0021 #include <cstdio>
0022 #include <map>
0023 #include <string>
0024 
0025 class TH2;
0026 class TH3;
0027 
0028 #define NumberOfOrders 15  // Convergence problems after 15; Rossegger used 30
0029 
0030 class Rossegger
0031 {
0032  public:
0033   explicit Rossegger(const std::string &filename);
0034   Rossegger(double a = 30, double b = 80, double L = 80, double epsilon = 1E-4);
0035   virtual ~Rossegger() {}
0036 
0037   void Verbosity(int v)
0038   {
0039     printf("verbosity set to %d.  was %d\n", v, verbosity);
0040     verbosity = v;
0041     return;
0042   };
0043   double Rmn(int m, int n, double r);               // Rmn function from Rossegger
0044   double Rmn_for_zeroes(int m, double x);           // Rmn function from Rossegger, as used to find Betamn zeroes.
0045   double Rmn1(int m, int n, double r);              // Rmn1 function from Rossegger
0046   double Rmn2(int m, int n, double r);              // Rmn2 function from Rossegger
0047   double RPrime(int m, int n, double a, double r);  // RPrime function from Rossegger
0048 
0049   double Rnk(int n, int k, double r);       // Rnk function from Rossegger
0050   double Rnk_for_zeroes(int n, double mu);  // Rnk function from Rossegger, as used to find munk zeroes.
0051 
0052   double Limu(double mu, double x);  // Bessel functions of purely imaginary order
0053   double Kimu(double mu, double x);  // Bessel functions of purely imaginary order
0054 
0055   double Ez(double r, double phi, double z, double r1, double phi1, double z1);
0056   double Er(double r, double phi, double z, double r1, double phi1, double z1);
0057   double Ephi(double r, double phi, double z, double r1, double phi1, double z1);
0058 
0059   // alternate versions that don't use precalc constants.
0060   double Rmn_(int m, int n, double r);  // Rmn function from Rossegger
0061   // Rmn_for_zeroes doesn't have a way to speed it up with precalcs.
0062   double Rmn1_(int m, int n, double r);              // Rmn1 function from Rossegger
0063   double Rmn2_(int m, int n, double r);              // Rmn2 function from Rossegger
0064   double RPrime_(int m, int n, double a, double r);  // RPrime function from Rossegger
0065   double Rnk_(int n, int k, double r);               // Rnk function from Rossegger
0066   double Rnk_for_zeroes_(int n, double mu);          // Rnk function from Rossegger, as used to find munk zeroes.
0067   double Ez_(double r, double phi, double z, double r1, double phi1, double z1);
0068   double Er_(double r, double phi, double z, double r1, double phi1, double z1);
0069   double Ephi_(double r, double phi, double z, double r1, double phi1, double z1);
0070 
0071  protected:
0072   bool fByFile = false;
0073   double a = NAN;
0074   double b = NAN;
0075   double L = NAN;  //  InnerRadius, OuterRadius, Length of 1/2 the TPC.
0076   int verbosity = 0;
0077   double pi = M_PI;
0078   double epsilon = NAN;  // precision.
0079 
0080   bool tweak = false;
0081 
0082   double MinimumDR = NAN;
0083   double MinimumDPHI = NAN;
0084   double MinimumDZ = NAN;
0085 
0086   double FindNextZero(double xstart, double epsilon, int order, double (Rossegger::*func)(int, double));  // Routine to find zeroes of func.
0087   void FindBetamn(double epsilon);                                                                        // Routine used to fill the Betamn array with resolution epsilon...
0088   void FindMunk(double epsilon);                                                                          // Routine used to fill the Munk array with resolution epsilon...
0089   bool CheckZeroes(double epsilon);                                                                       // confirm that the zeroes match to the desired precision.
0090 
0091   void LoadZeroes(const std::string &destfile);
0092   void SaveZeroes(const std::string &destfile);
0093 
0094   double Betamn[NumberOfOrders][NumberOfOrders]{};  //  Betamn array from Rossegger
0095   double N2mn[NumberOfOrders][NumberOfOrders]{};    //  N2mn array from Rossegger
0096   double Munk[NumberOfOrders][NumberOfOrders]{};    //  Munk array from Rossegger
0097   double N2nk[NumberOfOrders][NumberOfOrders]{};    //  N2nk array from Rossegger
0098 
0099   void PrecalcFreeConstants();  // Routine used to fill the arrays of other values used repeatedly:
0100   // helpful values to precompute for frequent use:
0101   double BetaN[NumberOfOrders]{};                               //  BetaN=(n+1)*pi/L as used in eg 5.32, .46
0102   double BetaN_a[NumberOfOrders]{};                             //  BetaN*a as used in eg 5.32, .46
0103   double BetaN_b[NumberOfOrders]{};                             //  BetaN*b as used in eg 5.32, .46
0104   double km_BetaN_a[NumberOfOrders][NumberOfOrders]{};          // kn(m,BetaN*a) as used in Rossegger 5.32
0105   double im_BetaN_a[NumberOfOrders][NumberOfOrders]{};          // in(m,BetaN*a) as used in Rossegger 5.32
0106   double km_BetaN_b[NumberOfOrders][NumberOfOrders]{};          // kn(m,BetaN*b) as used in Rossegger 5.33
0107   double im_BetaN_b[NumberOfOrders][NumberOfOrders]{};          // in(m,BetaN*b) as used in Rossegger 5.33
0108   double bessel_denominator[NumberOfOrders][NumberOfOrders]{};  // BesselI(m,BetaN*a)*BesselK(m,BetaN*b)-BesselI(m,BetaN*b)*BesselK(m,BetaN*a) as in Rossegger 5.65
0109 
0110   void PrecalcDerivedConstants();                           // Routine used to fill repeated values that depend on the Betamn and Munk zeroes:
0111   double ym_Betamn_a[NumberOfOrders][NumberOfOrders]{};     // yn(m,Betamn[m][n]*a) as used in Rossegger 5.11
0112   double jm_Betamn_a[NumberOfOrders][NumberOfOrders]{};     // jn(m,Betamn[m][n]*a) as used in Rossegger 5.11
0113   double liMunk_BetaN_a[NumberOfOrders][NumberOfOrders]{};  // limu(Munk[n][k],BetaN*a) as used in Rossegger 5.45
0114   double kiMunk_BetaN_a[NumberOfOrders][NumberOfOrders]{};  // kimu(Munk[n][k],BetaN*a) as used in Rossegger 5.45
0115   double sinh_Betamn_L[NumberOfOrders][NumberOfOrders]{};   // sinh(Betamn[m][n]*L)  as in Rossegger 5.64
0116   double sinh_pi_Munk[NumberOfOrders][NumberOfOrders]{};    // sinh(pi*Munk[n][k]) as in Rossegger 5.66
0117 
0118   TH2 *Tags = nullptr;
0119   std::map<std::string, TH3 *> Grid;
0120 };
0121 
0122 #endif /* __SPACECHARGE_H__ */