[82db145] | 1 | #include "TrkUtil.h"
|
---|
| 2 | #include <iostream>
|
---|
[c5696dd] | 3 | #include <algorithm>
|
---|
[82db145] | 4 |
|
---|
| 5 | // Constructor
|
---|
| 6 | TrkUtil::TrkUtil(Double_t Bz)
|
---|
| 7 | {
|
---|
| 8 | fBz = Bz;
|
---|
| 9 | fGasSel = 0; // Default is He-Isobuthane (90-10)
|
---|
| 10 | fRmin = 0.0; // Lower DCH radius
|
---|
| 11 | fRmax = 0.0; // Higher DCH radius
|
---|
| 12 | fZmin = 0.0; // Lower DCH z
|
---|
| 13 | fZmax = 0.0; // Higher DCH z
|
---|
| 14 | }
|
---|
| 15 | TrkUtil::TrkUtil()
|
---|
| 16 | {
|
---|
| 17 | fBz = 0.0;
|
---|
| 18 | fGasSel = 0; // Default is He-Isobuthane (90-10)
|
---|
| 19 | fRmin = 0.0; // Lower DCH radius
|
---|
| 20 | fRmax = 0.0; // Higher DCH radius
|
---|
| 21 | fZmin = 0.0; // Lower DCH z
|
---|
| 22 | fZmax = 0.0; // Higher DCH z
|
---|
| 23 | }
|
---|
| 24 | //
|
---|
| 25 | // Destructor
|
---|
| 26 | TrkUtil::~TrkUtil()
|
---|
| 27 | {
|
---|
| 28 | fBz = 0.0;
|
---|
| 29 | fGasSel = 0; // Default is He-Isobuthane (90-10)
|
---|
| 30 | fRmin = 0.0; // Lower DCH radius
|
---|
| 31 | fRmax = 0.0; // Higher DCH radius
|
---|
| 32 | fZmin = 0.0; // Lower DCH z
|
---|
| 33 | fZmax = 0.0; // Higher DCH z
|
---|
| 34 | }
|
---|
| 35 | //
|
---|
| 36 | // Helix parameters from position and momentum
|
---|
| 37 | // static
|
---|
| 38 | TVectorD TrkUtil::XPtoPar(TVector3 x, TVector3 p, Double_t Q, Double_t Bz)
|
---|
| 39 | {
|
---|
| 40 | //
|
---|
| 41 | TVectorD Par(5);
|
---|
| 42 | // Transverse parameters
|
---|
| 43 | Double_t a = -Q * Bz * cSpeed(); // Units are Tesla, GeV and meters
|
---|
| 44 | Double_t pt = p.Pt();
|
---|
| 45 | Double_t C = a / (2 * pt); // Half curvature
|
---|
| 46 | //std::cout << "ObsTrk::XPtoPar: fB = " << fB << ", a = " << a << ", pt = " << pt << ", C = " << C << std::endl;
|
---|
| 47 | Double_t r2 = x.Perp2();
|
---|
| 48 | Double_t cross = x(0) * p(1) - x(1) * p(0);
|
---|
[c5696dd] | 49 | Double_t T = sqrt(pt * pt - 2 * a * cross + a * a * r2);
|
---|
| 50 | Double_t phi0 = atan2((p(1) - a * x(0)) / T, (p(0) + a * x(1)) / T); // Phi0
|
---|
[82db145] | 51 | Double_t D; // Impact parameter D
|
---|
| 52 | if (pt < 10.0) D = (T - pt) / a;
|
---|
| 53 | else D = (-2 * cross + a * r2) / (T + pt);
|
---|
| 54 | //
|
---|
| 55 | Par(0) = D; // Store D
|
---|
| 56 | Par(1) = phi0; // Store phi0
|
---|
| 57 | Par(2) = C; // Store C
|
---|
| 58 | //Longitudinal parameters
|
---|
[c5696dd] | 59 | Double_t B = C * sqrt(TMath::Max(r2 - D * D, 0.0) / (1 + 2 * C * D));
|
---|
| 60 | Double_t st = asin(B) / C;
|
---|
[82db145] | 61 | Double_t ct = p(2) / pt;
|
---|
| 62 | Double_t z0 = x(2) - ct * st;
|
---|
| 63 | //
|
---|
| 64 | Par(3) = z0; // Store z0
|
---|
| 65 | Par(4) = ct; // Store cot(theta)
|
---|
| 66 | //
|
---|
| 67 | return Par;
|
---|
| 68 | }
|
---|
| 69 | // non-static
|
---|
| 70 | TVectorD TrkUtil::XPtoPar(TVector3 x, TVector3 p, Double_t Q)
|
---|
| 71 | {
|
---|
| 72 | //
|
---|
| 73 | TVectorD Par(5);
|
---|
| 74 | Double_t Bz = fBz;
|
---|
| 75 | Par = XPtoPar(x, p, Q, Bz);
|
---|
| 76 | //
|
---|
| 77 | return Par;
|
---|
| 78 | }
|
---|
| 79 | //
|
---|
| 80 | TVector3 TrkUtil::ParToX(TVectorD Par)
|
---|
| 81 | {
|
---|
| 82 | Double_t D = Par(0);
|
---|
| 83 | Double_t phi0 = Par(1);
|
---|
| 84 | Double_t z0 = Par(3);
|
---|
| 85 | //
|
---|
| 86 | TVector3 Xval;
|
---|
[c5696dd] | 87 | Xval(0) = -D * sin(phi0);
|
---|
| 88 | Xval(1) = D * cos(phi0);
|
---|
[82db145] | 89 | Xval(2) = z0;
|
---|
| 90 | //
|
---|
| 91 | return Xval;
|
---|
| 92 | }
|
---|
| 93 | //
|
---|
| 94 | TVector3 TrkUtil::ParToP(TVectorD Par)
|
---|
| 95 | {
|
---|
| 96 | if (fBz == 0.0)std::cout << "TrkUtil::ParToP: Warning Bz not set" << std::endl;
|
---|
| 97 | //
|
---|
| 98 | return ParToP(Par,fBz);
|
---|
| 99 | }
|
---|
| 100 | //
|
---|
| 101 | TVector3 TrkUtil::ParToP(TVectorD Par, Double_t Bz)
|
---|
| 102 | {
|
---|
| 103 | Double_t C = Par(2);
|
---|
| 104 | Double_t phi0 = Par(1);
|
---|
| 105 | Double_t ct = Par(4);
|
---|
| 106 | //
|
---|
| 107 | TVector3 Pval;
|
---|
| 108 | Double_t pt = Bz * cSpeed() / TMath::Abs(2 * C);
|
---|
[c5696dd] | 109 | Pval(0) = pt * cos(phi0);
|
---|
| 110 | Pval(1) = pt * sin(phi0);
|
---|
[82db145] | 111 | Pval(2) = pt * ct;
|
---|
| 112 | //
|
---|
| 113 | return Pval;
|
---|
| 114 | }
|
---|
| 115 | //
|
---|
| 116 | Double_t TrkUtil::ParToQ(TVectorD Par)
|
---|
| 117 | {
|
---|
| 118 | return TMath::Sign(1.0, -Par(2));
|
---|
| 119 | }
|
---|
| 120 |
|
---|
| 121 | //
|
---|
| 122 | // Parameter conversion to ACTS format
|
---|
| 123 | TVectorD TrkUtil::ParToACTS(TVectorD Par)
|
---|
| 124 | {
|
---|
| 125 | TVectorD pACTS(6); // Return vector
|
---|
| 126 | //
|
---|
| 127 | Double_t b = -cSpeed() * fBz / 2.;
|
---|
| 128 | pACTS(0) = 1000 * Par(0); // D from m to mm
|
---|
| 129 | pACTS(1) = 1000 * Par(3); // z0 from m to mm
|
---|
| 130 | pACTS(2) = Par(1); // Phi0 is unchanged
|
---|
[c5696dd] | 131 | pACTS(3) = atan2(1.0, Par(4)); // Theta in [0, pi] range
|
---|
| 132 | pACTS(4) = Par(2) / (b * sqrt(1 + Par(4) * Par(4))); // q/p in GeV
|
---|
[82db145] | 133 | pACTS(5) = 0.0; // Time: currently undefined
|
---|
| 134 | //
|
---|
| 135 | return pACTS;
|
---|
| 136 | }
|
---|
| 137 | // Covariance conversion to ACTS format
|
---|
| 138 | TMatrixDSym TrkUtil::CovToACTS(TVectorD Par, TMatrixDSym Cov)
|
---|
| 139 | {
|
---|
| 140 | TMatrixDSym cACTS(6); cACTS.Zero();
|
---|
| 141 | Double_t b = -cSpeed() * fBz / 2.;
|
---|
| 142 | //
|
---|
| 143 | // Fill derivative matrix
|
---|
| 144 | TMatrixD A(5, 5); A.Zero();
|
---|
| 145 | Double_t ct = Par(4); // cot(theta)
|
---|
| 146 | Double_t C = Par(2); // half curvature
|
---|
| 147 | A(0, 0) = 1000.; // D-D conversion to mm
|
---|
| 148 | A(1, 2) = 1.0; // phi0-phi0
|
---|
[c5696dd] | 149 | A(2, 4) = 1.0 / (sqrt(1.0 + ct * ct) * b); // q/p-C
|
---|
[82db145] | 150 | A(3, 1) = 1000.; // z0-z0 conversion to mm
|
---|
| 151 | A(4, 3) = -1.0 / (1.0 + ct * ct); // theta - cot(theta)
|
---|
| 152 | A(4, 4) = -C * ct / (b * pow(1.0 + ct * ct, 3.0 / 2.0)); // q/p-cot(theta)
|
---|
| 153 | //
|
---|
| 154 | TMatrixDSym Cv = Cov;
|
---|
| 155 | TMatrixD At(5, 5);
|
---|
| 156 | At.Transpose(A);
|
---|
| 157 | Cv.Similarity(At);
|
---|
| 158 | TMatrixDSub(cACTS, 0, 4, 0, 4) = Cv;
|
---|
| 159 | cACTS(5, 5) = 0.1; // Currently undefined: set to arbitrary value to avoid crashes
|
---|
| 160 | //
|
---|
| 161 | return cACTS;
|
---|
| 162 | }
|
---|
| 163 | //
|
---|
| 164 | // Parameter conversion to ILC format
|
---|
| 165 | TVectorD TrkUtil::ParToILC(TVectorD Par)
|
---|
| 166 | {
|
---|
| 167 | TVectorD pILC(5); // Return vector
|
---|
| 168 | //
|
---|
| 169 | pILC(0) = Par(0) * 1.0e3; // d0 in mm
|
---|
| 170 | pILC(1) = Par(1); // phi0 is unchanged
|
---|
| 171 | pILC(2) = -2 * Par(2) * 1.0e-3; // w in mm^-1
|
---|
| 172 | pILC(3) = Par(3) * 1.0e3; // z0 in mm
|
---|
| 173 | pILC(4) = Par(4); // tan(lambda) = cot(theta)
|
---|
| 174 | //
|
---|
| 175 | return pILC;
|
---|
| 176 | }
|
---|
| 177 | // Covariance conversion to ILC format
|
---|
| 178 | TMatrixDSym TrkUtil::CovToILC(TMatrixDSym Cov)
|
---|
| 179 | {
|
---|
| 180 | TMatrixDSym cILC(5); cILC.Zero();
|
---|
| 181 | //
|
---|
| 182 | // Fill derivative matrix
|
---|
| 183 | TMatrixD A(5, 5); A.Zero();
|
---|
| 184 | //
|
---|
| 185 | A(0, 0) = 1.0e3; // D-d0 in mm
|
---|
| 186 | A(1, 1) = 1.0; // phi0-phi0
|
---|
| 187 | A(2, 2) = -2.0e-3; // w-C
|
---|
| 188 | A(3, 3) = 1.0e3; // z0-z0 conversion to mm
|
---|
| 189 | A(4, 4) = 1.0; // tan(lambda) - cot(theta)
|
---|
| 190 | //
|
---|
| 191 | TMatrixDSym Cv = Cov;
|
---|
| 192 | TMatrixD At(5, 5);
|
---|
| 193 | At.Transpose(A);
|
---|
| 194 | Cv.Similarity(At);
|
---|
| 195 | cILC = Cv;
|
---|
| 196 | //
|
---|
| 197 | return cILC;
|
---|
| 198 | }
|
---|
| 199 | //
|
---|
| 200 | // Conversion from meters to mm
|
---|
| 201 | TVectorD TrkUtil::ParToMm(TVectorD Par) // Parameter conversion
|
---|
| 202 | {
|
---|
| 203 | TVectorD Pmm(5); // Return vector
|
---|
| 204 | //
|
---|
| 205 | Pmm(0) = Par(0) * 1.0e3; // d0 in mm
|
---|
| 206 | Pmm(1) = Par(1); // phi0 is unchanged
|
---|
| 207 | Pmm(2) = Par(2) * 1.0e-3; // C in mm^-1
|
---|
| 208 | Pmm(3) = Par(3) * 1.0e3; // z0 in mm
|
---|
| 209 | Pmm(4) = Par(4); // tan(lambda) = cot(theta) unchanged
|
---|
| 210 | //
|
---|
| 211 | return Pmm;
|
---|
| 212 | }
|
---|
| 213 | TMatrixDSym TrkUtil::CovToMm(TMatrixDSym Cov) // Covariance conversion
|
---|
| 214 | {
|
---|
| 215 | TMatrixDSym Cmm(5); Cmm.Zero();
|
---|
| 216 | //
|
---|
| 217 | // Fill derivative matrix
|
---|
| 218 | TMatrixD A(5, 5); A.Zero();
|
---|
| 219 | //
|
---|
| 220 | A(0, 0) = 1.0e3; // D-d0 in mm
|
---|
| 221 | A(1, 1) = 1.0; // phi0-phi0
|
---|
| 222 | A(2, 2) = 1.0e-3; // C-C
|
---|
| 223 | A(3, 3) = 1.0e3; // z0-z0 conversion to mm
|
---|
| 224 | A(4, 4) = 1.0; // lambda - cot(theta)
|
---|
| 225 | //
|
---|
| 226 | TMatrixDSym Cv = Cov;
|
---|
| 227 | TMatrixD At(5, 5);
|
---|
| 228 | At.Transpose(A);
|
---|
| 229 | Cv.Similarity(At);
|
---|
| 230 | Cmm = Cv;
|
---|
| 231 | //
|
---|
| 232 | return Cmm;
|
---|
| 233 | }
|
---|
| 234 | //
|
---|
[a95da74] | 235 |
|
---|
| 236 | //
|
---|
| 237 | void TrkUtil::SetBfield(Double_t Bz)
|
---|
| 238 | {
|
---|
| 239 | fBz = Bz;
|
---|
| 240 | }
|
---|
| 241 |
|
---|
[82db145] | 242 | // Setup chamber volume
|
---|
| 243 | void TrkUtil::SetDchBoundaries(Double_t Rmin, Double_t Rmax, Double_t Zmin, Double_t Zmax)
|
---|
| 244 | {
|
---|
| 245 | fRmin = Rmin; // Lower DCH radius
|
---|
| 246 | fRmax = Rmax; // Higher DCH radius
|
---|
| 247 | fZmin = Zmin; // Lower DCH z
|
---|
| 248 | fZmax = Zmax; // Higher DCH z
|
---|
| 249 | }
|
---|
| 250 | //
|
---|
| 251 | // Get Trakck length inside DCH volume
|
---|
| 252 | Double_t TrkUtil::TrkLen(TVectorD Par)
|
---|
| 253 | {
|
---|
| 254 | Double_t tLength = 0.0;
|
---|
| 255 | // Check if geometry is initialized
|
---|
| 256 | if (fZmin == 0.0 && fZmax == 0.0)
|
---|
| 257 | {
|
---|
| 258 | // No geometry set so send a warning and return 0
|
---|
| 259 | std::cout << "TrkUtil::TrkLen() called without a DCH volume defined" << std::endl;
|
---|
| 260 | }
|
---|
| 261 | else
|
---|
| 262 | {
|
---|
| 263 | //******************************************************************
|
---|
| 264 | // Determine the track length inside the chamber ****
|
---|
| 265 | //******************************************************************
|
---|
| 266 | //
|
---|
| 267 | // Track pararameters
|
---|
| 268 | Double_t D = Par(0); // Transverse impact parameter
|
---|
| 269 | Double_t phi0 = Par(1); // Transverse direction at minimum approach
|
---|
| 270 | Double_t C = Par(2); // Half curvature
|
---|
| 271 | Double_t z0 = Par(3); // Z at minimum approach
|
---|
| 272 | Double_t ct = Par(4); // cot(theta)
|
---|
| 273 | //std::cout << "TrkUtil:: parameters: D= " << D << ", phi0= " << phi0
|
---|
| 274 | // << ", C= " << C << ", z0= " << z0 << ", ct= " << ct << std::endl;
|
---|
| 275 | //
|
---|
[a95da74] | 276 | // Track length per unit phase change
|
---|
[c5696dd] | 277 | Double_t Scale = sqrt(1.0 + ct*ct) / (2.0*TMath::Abs(C));
|
---|
[82db145] | 278 | //
|
---|
| 279 | // Find intersections with chamber boundaries
|
---|
| 280 | //
|
---|
[a95da74] | 281 | Double_t phRin = 0.0; // phase of inner cylinder
|
---|
[82db145] | 282 | Double_t phRin2= 0.0; // phase of inner cylinder intersection (2nd branch)
|
---|
| 283 | Double_t phRhi = 0.0; // phase of outer cylinder intersection
|
---|
| 284 | Double_t phZmn = 0.0; // phase of left wall intersection
|
---|
| 285 | Double_t phZmx = 0.0; // phase of right wall intersection
|
---|
| 286 | // ... with inner cylinder
|
---|
| 287 | Double_t Rtop = TMath::Abs((1.0 + C*D) / C);
|
---|
| 288 |
|
---|
| 289 | if (Rtop > fRmin && TMath::Abs(D) < fRmin) // *** don't treat large D tracks for the moment ***
|
---|
| 290 | {
|
---|
[c5696dd] | 291 | Double_t ph = 2 * asin(C*sqrt((fRmin*fRmin - D*D) / (1.0 + 2.0*C*D)));
|
---|
[82db145] | 292 | Double_t z = z0 + ct*ph / (2.0*C);
|
---|
| 293 |
|
---|
| 294 | //std::cout << "Rin intersection: ph = " << ph<<", z= "<<z << std::endl;
|
---|
| 295 |
|
---|
[a95da74] | 296 | if (z < fZmax && z > fZmin) phRin = TMath::Abs(ph); // Intersection inside chamber volume
|
---|
[82db145] | 297 | //
|
---|
| 298 | // Include second branch of loopers
|
---|
[59ba063] | 299 | Double_t Pi = 3.14159265358979323846;
|
---|
| 300 | Double_t ph2 = 2*Pi - TMath::Abs(ph);
|
---|
[82db145] | 301 | if (ph < 0)ph2 = -ph2;
|
---|
| 302 | z = z0 + ct * ph2 / (2.0 * C);
|
---|
| 303 | if (z < fZmax && z > fZmin) phRin2 = TMath::Abs(ph2); // Intersection inside chamber volume
|
---|
| 304 | }
|
---|
| 305 | // ... with outer cylinder
|
---|
| 306 | if (Rtop > fRmax && TMath::Abs(D) < fRmax) // *** don't treat large D tracks for the moment ***
|
---|
| 307 | {
|
---|
[c5696dd] | 308 | Double_t ph = 2 * asin(C*sqrt((fRmax*fRmax - D*D) / (1.0 + 2.0*C*D)));
|
---|
[82db145] | 309 | Double_t z = z0 + ct*ph / (2.0*C);
|
---|
[a95da74] | 310 | if (z < fZmax && z > fZmin) phRhi = TMath::Abs(ph); // Intersection inside chamber volume
|
---|
[82db145] | 311 | }
|
---|
| 312 | // ... with left wall
|
---|
| 313 | Double_t Zdir = (fZmin - z0) / ct;
|
---|
| 314 | if (Zdir > 0.0)
|
---|
| 315 | {
|
---|
| 316 | Double_t ph = 2.0*C*Zdir;
|
---|
[c5696dd] | 317 | Double_t Rint = sqrt(D*D + (1.0 + 2.0*C*D)*pow(sin(ph / 2), 2) / (C*C));
|
---|
[a95da74] | 318 | if (Rint < fRmax && Rint > fRmin) phZmn = TMath::Abs(ph); // Intersection inside chamber volume
|
---|
[82db145] | 319 | }
|
---|
| 320 | // ... with right wall
|
---|
| 321 | Zdir = (fZmax - z0) / ct;
|
---|
| 322 | if (Zdir > 0.0)
|
---|
| 323 | {
|
---|
| 324 | Double_t ph = 2.0*C*Zdir;
|
---|
[c5696dd] | 325 | Double_t Rint = sqrt(D*D + (1.0 + 2.0*C*D)*pow(sin(ph / 2), 2) / (C*C));
|
---|
[a95da74] | 326 | if (Rint < fRmax && Rint > fRmin) phZmx = TMath::Abs(ph); // Intersection inside chamber volume
|
---|
[82db145] | 327 | }
|
---|
| 328 | //
|
---|
| 329 | // Order phases and keep the lowest two non-zero ones
|
---|
| 330 | //
|
---|
| 331 | const Int_t Nint = 5;
|
---|
| 332 | Double_t dPhase = 0.0; // Phase difference between two close intersections
|
---|
| 333 | Double_t ph_arr[Nint] = { phRin, phRin2, phRhi, phZmn, phZmx };
|
---|
[59ba063] | 334 | std::sort(ph_arr, ph_arr + Nint);
|
---|
[82db145] | 335 | Int_t iPos = -1; // First element > 0
|
---|
| 336 | for (Int_t i = 0; i < Nint; i++)
|
---|
| 337 | {
|
---|
[59ba063] | 338 | if (ph_arr[i] <= 0.0) iPos = i;
|
---|
[82db145] | 339 | }
|
---|
| 340 |
|
---|
| 341 | if (iPos < Nint - 2)
|
---|
| 342 | {
|
---|
[59ba063] | 343 | dPhase = ph_arr[iPos + 2] - ph_arr[iPos + 1];
|
---|
[82db145] | 344 | tLength = dPhase*Scale;
|
---|
| 345 | }
|
---|
| 346 | }
|
---|
| 347 | return tLength;
|
---|
| 348 | }
|
---|
| 349 | //
|
---|
| 350 | // Return number of ionization clusters
|
---|
| 351 | Bool_t TrkUtil::IonClusters(Double_t &Ncl, Double_t mass, TVectorD Par)
|
---|
| 352 | {
|
---|
| 353 | //
|
---|
| 354 | // Units are meters/Tesla/GeV
|
---|
| 355 | //
|
---|
| 356 | Ncl = 0.0;
|
---|
| 357 | Bool_t Signal = kFALSE;
|
---|
| 358 | Double_t tLen = 0;
|
---|
| 359 | // Check if geometry is initialized
|
---|
| 360 | if (fZmin == 0.0 && fZmax == 0.0)
|
---|
| 361 | {
|
---|
| 362 | // No geometry set so send a warning and return 0
|
---|
| 363 | std::cout << "TrkUtil::IonClusters() called without a volume defined" << std::endl;
|
---|
| 364 | }
|
---|
| 365 | else tLen = TrkLen(Par);
|
---|
| 366 |
|
---|
| 367 | //******************************************************************
|
---|
| 368 | // Now get the number of clusters ****
|
---|
| 369 | //******************************************************************
|
---|
| 370 | //
|
---|
| 371 | Double_t muClu = 0.0; // mean number of clusters
|
---|
| 372 | Double_t bg = 0.0; // beta*gamma
|
---|
| 373 | Ncl = 0.0;
|
---|
| 374 | if (tLen > 0.0)
|
---|
| 375 | {
|
---|
| 376 | Signal = kTRUE;
|
---|
| 377 | //
|
---|
| 378 | // Find beta*gamma
|
---|
| 379 | if (fBz == 0.0)
|
---|
| 380 | {
|
---|
| 381 | Signal = kFALSE;
|
---|
| 382 | std::cout << "TrkUtil::IonClusters: Please set Bz!!!" << std::endl;
|
---|
| 383 | }
|
---|
| 384 | else
|
---|
| 385 | {
|
---|
| 386 | TVector3 p = ParToP(Par);
|
---|
| 387 | bg = p.Mag() / mass;
|
---|
| 388 | muClu = Nclusters(bg)*tLen; // Avg. number of clusters
|
---|
| 389 |
|
---|
| 390 | Ncl = gRandom->PoissonD(muClu); // Actual number of clusters
|
---|
| 391 | }
|
---|
| 392 |
|
---|
| 393 | }
|
---|
| 394 | //
|
---|
| 395 | return Signal;
|
---|
| 396 | }
|
---|
| 397 | //
|
---|
| 398 | //
|
---|
[a95da74] | 399 | Double_t TrkUtil::Nclusters(Double_t begam)
|
---|
[82db145] | 400 | {
|
---|
| 401 | Int_t Opt = fGasSel;
|
---|
| 402 | Double_t Nclu = Nclusters(begam, Opt);
|
---|
| 403 | //
|
---|
| 404 | return Nclu;
|
---|
| 405 | }
|
---|
| 406 | //
|
---|
| 407 | Double_t TrkUtil::Nclusters(Double_t begam, Int_t Opt) {
|
---|
| 408 | //
|
---|
| 409 | // Opt = 0: He 90 - Isobutane 10
|
---|
| 410 | // = 1: pure He
|
---|
| 411 | // = 2: Argon 50 - Ethane 50
|
---|
| 412 | // = 3: pure Argon
|
---|
| 413 | //
|
---|
| 414 | //
|
---|
[4df491e] | 415 | /*
|
---|
[82db145] | 416 | std::vector<double> bg{ 0.5, 0.8, 1., 2., 3., 4., 5., 8., 10.,
|
---|
| 417 | 12., 15., 20., 50., 100., 200., 500., 1000. };
|
---|
| 418 | // He 90 - Isobutane 10
|
---|
| 419 | std::vector<double> ncl_He_Iso{ 42.94, 23.6,18.97,12.98,12.2,12.13,
|
---|
| 420 | 12.24,12.73,13.03,13.29,13.63,14.08,15.56,16.43,16.8,16.95,16.98 };
|
---|
| 421 | //
|
---|
| 422 | // pure He
|
---|
| 423 | std::vector<double> ncl_He{ 11.79,6.5,5.23,3.59,3.38,3.37,3.4,3.54,3.63,
|
---|
| 424 | 3.7,3.8,3.92,4.33,4.61,4.78,4.87,4.89 };
|
---|
| 425 | //
|
---|
| 426 | // Argon 50 - Ethane 50
|
---|
| 427 | std::vector<double> ncl_Ar_Eth{ 130.04,71.55,57.56,39.44,37.08,36.9,
|
---|
| 428 | 37.25,38.76,39.68,40.49,41.53,42.91,46.8,48.09,48.59,48.85,48.93 };
|
---|
| 429 | //
|
---|
| 430 | // pure Argon
|
---|
| 431 | std::vector<double> ncl_Ar{ 88.69,48.93,39.41,27.09,25.51,25.43,25.69,
|
---|
| 432 | 26.78,27.44,28.02,28.77,29.78,32.67,33.75,34.24,34.57,34.68 };
|
---|
| 433 | //
|
---|
| 434 | Int_t nPoints = (Int_t)bg.size();
|
---|
| 435 | bg.push_back(10000.);
|
---|
| 436 | std::vector<double> ncl;
|
---|
| 437 | switch (Opt)
|
---|
| 438 | {
|
---|
| 439 | case 0: ncl = ncl_He_Iso; // He-Isobutane
|
---|
| 440 | break;
|
---|
| 441 | case 1: ncl = ncl_He; // pure He
|
---|
| 442 | break;
|
---|
| 443 | case 2: ncl = ncl_Ar_Eth; // Argon - Ethane
|
---|
| 444 | break;
|
---|
| 445 | case 3: ncl = ncl_Ar; // pure Argon
|
---|
| 446 | break;
|
---|
| 447 | }
|
---|
| 448 | ncl.push_back(ncl[nPoints - 1]);
|
---|
[4df491e] | 449 | */
|
---|
| 450 | const Int_t Npt = 18;
|
---|
| 451 | Double_t bg[Npt] = { 0.5, 0.8, 1., 2., 3., 4., 5., 8., 10.,
|
---|
| 452 | 12., 15., 20., 50., 100., 200., 500., 1000., 10000. };
|
---|
| 453 | //
|
---|
| 454 | // He 90 - Isobutane 10
|
---|
| 455 | Double_t ncl_He_Iso[Npt] = { 42.94, 23.6,18.97,12.98,12.2,12.13,
|
---|
| 456 | 12.24,12.73,13.03,13.29,13.63,14.08,15.56,16.43,16.8,16.95,16.98, 16.98 };
|
---|
| 457 | //
|
---|
| 458 | // pure He
|
---|
| 459 | Double_t ncl_He[Npt] = { 11.79,6.5,5.23,3.59,3.38,3.37,3.4,3.54,3.63,
|
---|
| 460 | 3.7,3.8,3.92,4.33,4.61,4.78,4.87,4.89, 4.89 };
|
---|
| 461 | //
|
---|
| 462 | // Argon 50 - Ethane 50
|
---|
| 463 | Double_t ncl_Ar_Eth[Npt] = { 130.04,71.55,57.56,39.44,37.08,36.9,
|
---|
| 464 | 37.25,38.76,39.68,40.49,41.53,42.91,46.8,48.09,48.59,48.85,48.93,48.93 };
|
---|
| 465 | //
|
---|
| 466 | // pure Argon
|
---|
| 467 | Double_t ncl_Ar[Npt] = { 88.69,48.93,39.41,27.09,25.51,25.43,25.69,
|
---|
| 468 | 26.78,27.44,28.02,28.77,29.78,32.67,33.75,34.24,34.57,34.68, 34.68 };
|
---|
| 469 | //
|
---|
| 470 | Double_t ncl[Npt];
|
---|
| 471 | switch (Opt)
|
---|
| 472 | {
|
---|
| 473 | case 0: std::copy(ncl_He_Iso, ncl_He_Iso + Npt, ncl); // He-Isobutane
|
---|
[a95da74] | 474 | break;
|
---|
[4df491e] | 475 | case 1: std::copy(ncl_He, ncl_He + Npt, ncl); // pure He
|
---|
| 476 | break;
|
---|
| 477 | case 2: std::copy(ncl_Ar_Eth, ncl_Ar_Eth + Npt, ncl); // Argon - Ethane
|
---|
| 478 | break;
|
---|
| 479 | case 3: std::copy(ncl_Ar, ncl_Ar + Npt, ncl); // pure Argon
|
---|
| 480 | break;
|
---|
| 481 | }
|
---|
| 482 | //
|
---|
[82db145] | 483 | Int_t ilow = 0;
|
---|
| 484 | while (begam > bg[ilow])ilow++;
|
---|
| 485 | ilow--;
|
---|
| 486 | //std::cout << "ilow= " << ilow << ", low = " << bg[ilow] << ", val = " << begam
|
---|
| 487 | // << ", high = " << bg[ilow + 1] << std::endl;
|
---|
| 488 | //
|
---|
| 489 | Int_t ind[3] = { ilow, ilow + 1, ilow + 2 };
|
---|
| 490 | TVectorD y(3);
|
---|
| 491 | for (Int_t i = 0; i < 3; i++)y(i) = ncl[ind[i]];
|
---|
| 492 | TVectorD x(3);
|
---|
| 493 | for (Int_t i = 0; i < 3; i++)x(i) = bg[ind[i]];
|
---|
| 494 | TMatrixD Xval(3, 3);
|
---|
| 495 | for (Int_t i = 0; i < 3; i++)Xval(i, 0) = 1.0;
|
---|
| 496 | for (Int_t i = 0; i < 3; i++)Xval(i, 1) = x(i);
|
---|
| 497 | for (Int_t i = 0; i < 3; i++)Xval(i, 2) = x(i) * x(i);
|
---|
| 498 | //std::cout << "Xval:" << std::endl; Xval.Print();
|
---|
| 499 | Xval.Invert();
|
---|
| 500 | TVectorD coeff = Xval * y;
|
---|
| 501 | Double_t interp = coeff[0] + coeff[1] * begam + coeff[2] * begam * begam;
|
---|
[a95da74] | 502 | //std::cout << "val1= (" <<x(0)<<", "<< y(0) << "), val2= ("
|
---|
| 503 | // <<x(1)<<", "<< y(1) << "), val3= ("
|
---|
[82db145] | 504 | // <<x(2)<<", "<< y(2)
|
---|
| 505 | // << "), result= (" <<begam<<", "<< interp<<")" << std::endl;
|
---|
| 506 | //
|
---|
[c5696dd] | 507 | //if (TMath::IsNaN(interp))std::cout << "NaN found: bg= " << begam << ", Opt= " << Opt << std::endl;
|
---|
[82db145] | 508 | if (begam < bg[0]) interp = 0.0;
|
---|
| 509 | //std::cout << "bg= " << begam << ", Opt= " << Opt <<", interp = "<<interp<< std::endl;
|
---|
| 510 | return 100*interp;
|
---|
| 511 | }
|
---|
| 512 | //
|
---|
| 513 | Double_t TrkUtil::funcNcl(Double_t *xp, Double_t *par){
|
---|
| 514 | Double_t bg = xp[0];
|
---|
| 515 | return Nclusters(bg);
|
---|
| 516 | }
|
---|
| 517 | //
|
---|
| 518 | void TrkUtil::SetGasMix(Int_t Opt)
|
---|
| 519 | {
|
---|
| 520 | if (Opt < 0 || Opt > 3)
|
---|
| 521 | {
|
---|
| 522 | std::cout << "TrkUtil::SetGasMix Gas option not allowed. No action."
|
---|
| 523 | << std::endl;
|
---|
| 524 | }
|
---|
| 525 | else fGasSel = Opt;
|
---|
[c5696dd] | 526 | }
|
---|