Project
Loading...
Searching...
No Matches
mathieson.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
12#include <cstdio>
13#include <cstdlib>
14#include <stdexcept>
15#include <map>
16#include <limits>
17
19#include "mathUtil.h"
20#include "mathieson.h"
21
22namespace o2
23{
24namespace mch
25{
26extern ClusterConfig clusterConfig;
27
28// Chamber 1, 2
29const double sqrtK3x1_2 = 0.7000; // Pitch= 0.21 cm
30const double sqrtK3y1_2 = 0.7550; // Pitch= 0.21 cm
31const double pitch1_2 = 0.21;
32// Chamber 3, 10
33const double sqrtK3x3_10 = 0.7131; // Pitch= 0.25 cm
34const double sqrtK3y3_10 = 0.7642; // Pitch= 0.25 cm
35const double pitch3_10 = 0.25;
36
37int mathiesonType; // 0 for Station 1 or 1 for station 2-5
38static double K1x[2], K1y[2];
39static double K2x[2], K2y[2];
40static const double sqrtK3x[2] = {sqrtK3x1_2, sqrtK3x3_10},
41 sqrtK3y[2] = {sqrtK3y1_2, sqrtK3y3_10};
42static double K3x[2], K3y[2];
43static double K4x[2], K4y[2];
44static double pitch[2] = {pitch1_2, pitch3_10};
45static double invPitch[2];
46
47// Spline Coef
48int useSpline = 0;
50static double splineXYStep = 1.0e-3;
51static double splineXYLimit = 3.0;
52static int nSplineSampling = 0;
53double* splineXY = nullptr;
54
55//
56int useCache = 0;
57
59{
60 a = new double[N];
61 b = new double[N];
62 c = new double[N];
63 d = new double[N];
64}
65
67{
68 delete[] a;
69 delete[] b;
70 delete[] c;
71 delete[] d;
72}
73void initMathieson(int useSpline_, int useCache_)
74{
75 useSpline = useSpline_;
76 useCache = useCache_;
77 //
78 for (int i = 0; i < 2; i++) {
79 K3x[i] = sqrtK3x[i] * sqrtK3x[i];
80 K3y[i] = sqrtK3y[i] * sqrtK3y[i];
81 K2x[i] = M_PI * 0.5 * (1.0 - sqrtK3x[i] * 0.5);
82 K2y[i] = M_PI * 0.5 * (1.0 - sqrtK3y[i] * 0.5);
83 K1x[i] = K2x[i] * sqrtK3x[i] * 0.25 / (atan(sqrtK3x[i]));
84 K1y[i] = K2y[i] * sqrtK3y[i] * 0.25 / (atan(sqrtK3y[i]));
85 K4x[i] = K1x[i] / K2x[i] / sqrtK3x[i];
86 K4y[i] = K1y[i] / K2y[i] / sqrtK3y[i];
87 invPitch[i] = 1.0 / pitch[i];
88 }
89 if (useSpline) {
91 }
92}
93
95{
96 // x/y Interval and positive x/y limit
97 double xyStep = splineXYStep;
98 double xyLimit = splineXYLimit;
99 // X/Y Sampling
100 nSplineSampling = int(xyLimit / xyStep) + 1;
101 int N = nSplineSampling;
102
103 splineXY = new double[N];
104 for (int i = 0; i < N; i++) {
105 splineXY[i] = xyStep * i;
106 }
107 double* xy = splineXY;
108
109 // Spline coef allocation for the 4 functions
110 splineCoef[0][0] = new SplineCoef(N);
111 splineCoef[0][1] = new SplineCoef(N);
112 splineCoef[1][0] = new SplineCoef(N);
113 splineCoef[1][1] = new SplineCoef(N);
114
115 // Compute the spline Coef. for the 4 Mathieson primitives
116 double mathPrimitive[N];
117 double rightDerivative(0.0), leftDerivative;
118 // X and Y primitives on chambers <= 2 (Mathieson Type = 0)
119 int mathiesonType = 0;
120 int axe = 0;
121 mathiesonPrimitive(xy, N, axe, 2, mathPrimitive);
122 leftDerivative = 2.0 * K4x[mathiesonType] * sqrtK3x[mathiesonType] * K2x[mathiesonType] * invPitch[mathiesonType];
123 computeSplineCoef(xy, xyStep, mathPrimitive, N, leftDerivative, rightDerivative, o2::mch::splineCoef[mathiesonType][axe]);
124 axe = 1;
125 mathiesonPrimitive(xy, N, axe, 2, mathPrimitive);
126 leftDerivative = 2.0 * K4y[mathiesonType] * sqrtK3y[mathiesonType] * K2y[mathiesonType] * invPitch[mathiesonType];
127 computeSplineCoef(xy, xyStep, mathPrimitive, N, leftDerivative, rightDerivative, splineCoef[mathiesonType][axe]);
128 mathiesonType = 1;
129 axe = 0;
130 mathiesonPrimitive(xy, N, axe, 3, mathPrimitive);
131 leftDerivative = 2.0 * K4x[mathiesonType] * sqrtK3x[mathiesonType] * K2x[mathiesonType] * invPitch[mathiesonType];
132 computeSplineCoef(xy, xyStep, mathPrimitive, N, leftDerivative, rightDerivative, splineCoef[mathiesonType][axe]);
133 axe = 1;
134 mathiesonPrimitive(xy, N, axe, 3, mathPrimitive);
135 leftDerivative = 2.0 * K4y[mathiesonType] * sqrtK3y[mathiesonType] * K2y[mathiesonType] * invPitch[mathiesonType];
136 computeSplineCoef(xy, xyStep, mathPrimitive, N, leftDerivative, rightDerivative, splineCoef[mathiesonType][axe]);
137}
138
139// Spline implementation of the book "Numerical Analysis" - 9th edition
140// Richard L Burden, J Douglas Faires
141// Section 3.5, p. 146
142// Restrictions : planed with a regular sampling (dx = cst)
143// spline(x) :[-inf, +inf] -> [-1/2, +1/2]
144// Error < 7.0 e-11 for 1001 sampling between [0, 3.0]
145void computeSplineCoef(const double* xy, double xyStep, const double* f, int N,
146 double leftDerivative, double rightDerivative, SplineCoef* splineCoef)
147{
148 double* a = splineCoef->a;
149 double* b = splineCoef->b;
150 double* c = splineCoef->c;
151 double* d = splineCoef->d;
152
153 // a coef : the sampled function
154 vectorCopy(f, N, a);
155
156 // Step 1
157 double h = xyStep;
158
159 // Step 2 & 3 : Compute alpha
160 double alpha[N];
161 alpha[0] = 3.0 / h * (f[1] - f[0]) - 3 * leftDerivative;
162 alpha[N - 1] = 3 * rightDerivative - 3.0 / h * (f[N - 1] - f[N - 2]);
163 for (int i = 1; i < N - 1; i++) {
164 // To keep the general case if h is not constant
165 alpha[i] = 3.0 / h * (f[i + 1] - f[i]) - 3.0 / h * (f[i] - f[i - 1]);
166 }
167
168 // Step 4 to 6 solve a tridiagonal linear system
169 //
170 // Step 4
171 double l[N], mu[N], z[N];
172 l[0] = 2.0 * h;
173 mu[0] = 0.5;
174 z[0] = alpha[0] / l[0];
175 //
176 // Step 5
177 for (int i = 1; i < N - 1; i++) {
178 l[i] = 2.0 * (xy[i + 1] - xy[i - 1]) - h * mu[i - 1];
179 mu[i] = h / l[i];
180 z[i] = (alpha[i] - h * z[i - 1]) / l[i];
181 }
182 //
183 // Step 6
184 l[N - 1] = h * (2.0 - mu[N - 2]);
185 z[N - 1] = (alpha[N - 1] - h * z[N - 2]) / l[N - 1];
186 c[N - 1] = z[N - 1];
187
188 // Step 7 : calculate cubic coefficients
189 for (int j = N - 2; j >= 0; j--) {
190 c[j] = z[j] - mu[j] * c[j + 1];
191 b[j] = (f[j + 1] - f[j]) / h - h / 3.0 * (c[j + 1] + 2 * c[j]);
192 d[j] = (c[j + 1] - c[j]) / (3 * h);
193 }
194}
195
196void splineMathiesonPrimitive(const double* x, int N, int axe, int chamberId, double* mPrimitive)
197{
198 int mathiesonType = (chamberId <= 2) ? 0 : 1;
199 double* a = splineCoef[mathiesonType][axe]->a;
200 double* b = splineCoef[mathiesonType][axe]->b;
201 double* c = splineCoef[mathiesonType][axe]->c;
202 double* d = splineCoef[mathiesonType][axe]->d;
203 double dx = splineXYStep;
204 // printf("dx=%f nSplineSampling=%d\n", dx, nSplineSampling);
205 double signX[N];
206 // x without sign
207 double uX[N];
208 for (int i = 0; i < N; i++) {
209 signX[i] = (x[i] >= 0) ? 1 : -1;
210 uX[i] = signX[i] * x[i];
211 /*
212 if( uX[i] > (2.0 * splineXYLimit)) {
213 // x >> 0, f(x) = 0.0
214 signX[i] = 0.0;
215 uX[i] = 0.5;
216 }
217 */
218 }
219
220 double cst = 1.0 / dx;
221 // Get indexes in the sample function
222 int idx;
223 double h;
224 for (int i = 0; i < N; i++) {
225 // int k = int( uX[i] * cst + dx*0.1 );
226 // if ( k < nSplineSampling) {
227 if (uX[i] < splineXYLimit) {
228 idx = int(uX[i] * cst + dx * 0.1);
229 h = uX[i] - idx * dx;
230 } else {
231 idx = nSplineSampling - 1;
232 h = 0;
233 }
234 mPrimitive[i] = signX[i] * (a[idx] + h * (b[idx] + h * (c[idx] + h * (d[idx]))));
235 // printf("x[i]=%f, signX[i]=%f uX[i]=%f idx=%d, h=%f, prim=%f, splineXYLimit=%f\n", x[i], signX[i], uX[i], idx, h, mPrimitive[i], splineXYLimit );
236 }
237 // print ("uX ", uX)
238 // print ("h ", h)
239 // print ("f(x0) ", a[idx])
240 // print ("df|dx0", h*( b[idx] + h*( c[idx] + h *(d[idx]))))
241 // print ("f, ", a[idx] + h*( b[idx] + h*( c[idx] + h *(d[idx]))))
242}
243
244// Return the Mathieson primitive at x or y
245void mathiesonPrimitive(const double* xy, int N,
246 int axe, int chamberId, double mPrimitive[])
247{
248 mathiesonType = (chamberId <= 2) ? 0 : 1;
249 //
250 // Select Mathieson coef.
251 double curK2xy = (axe == 0) ? K2x[mathiesonType] : K2y[mathiesonType];
252 double curSqrtK3xy = (axe == 0) ? sqrtK3x[mathiesonType] : sqrtK3y[mathiesonType];
253 double curInvPitch = invPitch[mathiesonType];
254 double cst2xy = curK2xy * curInvPitch;
255 double curK4xy = (axe == 0) ? K4x[mathiesonType] : K4y[mathiesonType];
256
257 for (int i = 0; i < N; i++) {
258 double u = curSqrtK3xy * tanh(cst2xy * xy[i]);
259 mPrimitive[i] = 2 * curK4xy * atan(u);
260 }
261}
262
263void compute1DMathieson(const double* xy, int N,
264 int axe, int chamberId, double mathieson[])
265{
266 // Returning array: Charge Integral on all the pads
267 //
268 mathiesonType = (chamberId <= 2) ? 0 : 1;
269
270 //
271 // Select Mathieson coef.
272
273 double curK1xy = (axe == 0) ? K1x[mathiesonType] : K1y[mathiesonType];
274 double curK2xy = (axe == 0) ? K2x[mathiesonType] : K2y[mathiesonType];
275 double curK3xy = (axe == 0) ? K3x[mathiesonType] : K3y[mathiesonType];
276 double curInvPitch = invPitch[mathiesonType];
277 double cst2xy = curK2xy * curInvPitch;
278
279 for (int i = 0; i < N; i++) {
280 // tanh(x) & tanh(y)
281 double xTanh = tanh(cst2xy * xy[i]);
282 double xTanh2 = xTanh * xTanh;
283 mathieson[i] = curK1xy * (1.0 - xTanh2) / (1.0 + curK3xy * xTanh2);
284 }
285 return;
286}
287void compute1DPadIntegrals(const double* xyInf, const double* xySup, int N,
288 double xy0, int axe, int chamberId, double* integrals)
289{
290 double zInf[N], zSup[N];
291 vectorAddScalar(xyInf, -xy0, N, zInf);
292 vectorAddScalar(xySup, -xy0, N, zSup);
293 compute1DPadIntegrals(zInf, zSup, N, axe, chamberId, integrals);
294}
295
296void compute1DPadIntegrals(const double* xyInf, const double* xySup, int N,
297 int axe, int chamberId, double* Integrals)
298{
299 // Returning array: Charge Integral on all the pads
300 //
301 mathiesonType = (chamberId <= 2) ? 0 : 1;
302
303 //
304 // Select Mathieson coef.
305 double curInvPitch = invPitch[mathiesonType];
306 double curK2 = (axe == 0) ? K2x[mathiesonType] : K2y[mathiesonType];
307 double curSqrtK3 = (axe == 0) ? sqrtK3x[mathiesonType] : sqrtK3y[mathiesonType];
308 double curK4 = (axe == 0) ? K4x[mathiesonType] : K4y[mathiesonType];
309 double cst2 = curK2 * curInvPitch;
310 double cst4 = 2.0 * curK4;
311
312 double uInf, uSup;
313 for (int i = 0; i < N; i++) {
314 // x/u
315 uInf = curSqrtK3 * tanh(cst2 * xyInf[i]);
316 uSup = curSqrtK3 * tanh(cst2 * xySup[i]);
317 //
318 Integrals[i] = cst4 * (atan(uSup) - atan(uInf));
319 // printf(" xyInfSup %2d [%10.6g, %10.6g] x [%10.6g, %10.6g]-> %10.6g\n",
320 // i, xInf[i], xSup[i], yInf[i], ySup[i], Integrals[i]);
321 }
322 // printf(" I[0..%3ld] = %f, %f, ... %f\n", N-1, Integrals[0], Integrals[1],
323 // Integrals[N-1]);
324 return;
325}
326
327int compressSameValues(const double* x1, const double* x2, int* map1, int* map2, int N, double* xCompress)
328{
329 // map1[0..N-1]: i in [0..N-1] -> integral index for x1 [0..nCompressed-1]
330 // map2[0..N-1]: the same for x2
331 // xCompress[0..nCompressed]: values of x1 & x2 compressed (unique values)
332 // The xCompress values will be used to compute the primitive
333 // The map1/2 will be used to find the corresponding index in the xCompress or primitive arrays
334 // Return nCompressed
335
336 // Transform to integer to avoid comparison on close x values
337 const double* x[2] = {x1, x2};
338 int* xCode = new int[2 * N];
339 for (int i = 0; i < N; i++) {
340 for (int b = 0; b < 2; b++) {
341 // Calculate the indexes in the 1D charge integral
342 // Error on pad position > 10-3 cm
343 xCode[i + b * N] = (int)(x[b][i] * 1000 + 0.5);
344 }
345 }
346 // Sort the code
347 int sIdx[2 * N];
348 for (int k = 0; k < 2 * N; k++) {
349 sIdx[k] = k;
350 }
351 std::sort(sIdx, &sIdx[2 * N], [=](int a, int b) {
352 return (xCode[a] < xCode[b]);
353 });
354
355 // printf("sort xCode[sIdx[0]]=%d xCode[sIdx[2*N-1]]=%d\n", xCode[sIdx[0]], xCode[sIdx[2*N-1]]);
356 // vectorPrintInt("xCode",xCode, 2*N);
357 // vectorPrintInt("sIdx",sIdx, 2*N);
358
359 // Renumber and compress
360 int nCompress = 0;
361 int prevCode = std::numeric_limits<int>::max();
362
363 // Map1
364 for (int i = 0; i < 2 * N; i++) {
365 int idx = sIdx[i];
366 if (xCode[idx] != prevCode) {
367 if (idx < N) {
368 // Store the compress value in map1
369 xCompress[nCompress] = x1[idx];
370 map1[idx] = nCompress;
371 // printf("i=%d sIdx[i]=%d nCompress=%d idx=%d map1[idx]=%d\n", i, idx, nCompress, idx, map1[idx]);
372 } else {
373 // Store the compress value in map2
374 xCompress[nCompress] = x2[idx - N];
375 map2[idx - N] = nCompress;
376 // printf("i=%d sIdx[i]=%d nCompress=%d idx-N=%d map2[idx]=%d\n", i, idx, nCompress, idx-N, map2[idx-N]);
377 }
378 nCompress++;
379 } else {
380 // the code is the same (same values)
381 if (idx < N) {
382 map1[idx] = nCompress - 1;
383 // printf("identical i=%d sIdx[i]=%d nCompress-1=%d idx=%d\n", i, idx, nCompress-1, idx);
384 } else {
385 map2[idx - N] = nCompress - 1;
386 // printf("identical i=%d sIdx[i]=%d nCompress-1=%d idx=%d\n", i, idx, nCompress-1, idx-N);
387 }
388 }
389 prevCode = xCode[idx];
390 }
391 // printf(" compress nCompress/N=%d/%d \n", nCompress, N);
392 // vectorPrint("x1", x1, N);
393 // vectorPrintInt("map1",map1, N);
394 // vectorPrint("x2", x2, N);
395 // vectorPrintInt("map2",map2, N);
396 // vectorPrint("xCompress", xCompress, nCompress);
397 delete[] xCode;
398 return nCompress;
399}
400
401CompressedPads_t* compressPads(const double* xInf, const double* xSup,
402 const double* yInf, const double* ySup, int N)
403{
404 CompressedPads_t* compressedPads = new CompressedPads_t;
405 // On x axe
406 compressedPads->xCompressed = new double[2 * N];
407 compressedPads->mapXInf = new int[N];
408 compressedPads->mapXSup = new int[N];
409 compressedPads->nXc = compressSameValues(xInf, xSup, compressedPads->mapXInf, compressedPads->mapXSup, N, compressedPads->xCompressed);
410 compressedPads->yCompressed = new double[2 * N];
411 compressedPads->mapYInf = new int[N];
412 compressedPads->mapYSup = new int[N];
413 compressedPads->nYc = compressSameValues(yInf, ySup, compressedPads->mapYInf, compressedPads->mapYSup, N, compressedPads->yCompressed);
414 return compressedPads;
415}
416
418{
419 delete[] compressedPads->mapXInf;
420 delete[] compressedPads->mapXSup;
421 delete[] compressedPads->mapYInf;
422 delete[] compressedPads->mapYSup;
423 delete[] compressedPads->xCompressed;
424 delete[] compressedPads->yCompressed;
425}
426
428 /* const double* xInf, const double* xSup,
429 const double* yInf, const double* ySup,
430 */
431 CompressedPads_t* compressedPads, double xShift, double yShift, int N,
432 int chamberId, double Integrals[])
433{
434
435 int nXc = compressedPads->nXc;
436 int nYc = compressedPads->nYc;
437 // Compute the integrals on Compressed pads
438 double xy[N];
439 double xPrimitives[nXc];
440 double yPrimitives[nYc];
441 // X axe
442 int axe = 0;
443 // x Translation (seed location)
444 vectorAddScalar(compressedPads->xCompressed, -xShift, nXc, xy);
445 // Primitives on compressed pads
446 mathiesonPrimitive(xy, nXc, axe, chamberId, xPrimitives);
447 // Y axe
448 axe = 1;
449 // x Translation (seed location)
450 vectorAddScalar(compressedPads->yCompressed, -yShift, nYc, xy);
451 // Primitives on compressed pads
452 mathiesonPrimitive(xy, nYc, axe, chamberId, yPrimitives);
453
454 // Compute all the integrals
455 int* mapXInf = compressedPads->mapXInf;
456 int* mapXSup = compressedPads->mapXSup;
457 int* mapYInf = compressedPads->mapYInf;
458 int* mapYSup = compressedPads->mapYSup;
459 for (int i = 0; i < N; i++) {
460 Integrals[i] = (xPrimitives[mapXSup[i]] - xPrimitives[mapXInf[i]]) * (yPrimitives[mapYSup[i]] - yPrimitives[mapYInf[i]]);
461 // printf(" i=%d mapXInf=%d mapXSup=%d mapYInf=%d mapYSup=%d xyIntegrals=%f, %f \n", i,
462 // mapXInf[i], mapXSup[i], mapYInf[i], mapYSup[i], xPrimitives[mapXSup[i]] - xPrimitives[mapXInf[i]],
463 // yPrimitives[mapYSup[i]] - yPrimitives[mapYInf[i]]);
464 }
465
466 // vectorPrint("xPrimitives", xPrimitives, nXc);
467 // vectorPrint("yPrimitives", yPrimitives, nYc);
468}
469
470void compute2DPadIntegrals(const double* xInf, const double* xSup,
471 const double* yInf, const double* ySup, int N,
472 int chamberId, double Integrals[])
473{
474 if (1) {
475 int mapXInf[N], mapXSup[N];
476 int mapYInf[N], mapYSup[N];
477 double xy[2 * N];
478 // Primitives on x axe
479 int nXc = compressSameValues(xInf, xSup, mapXInf, mapXSup, N, xy);
480 // vectorPrint("x map", xy, nXc);
481 int axe = 0;
482 double xPrimitives[nXc];
483 mathiesonPrimitive(xy, nXc, axe, chamberId, xPrimitives);
484 // Primitives on y axe
485 int nYc = compressSameValues(yInf, ySup, mapYInf, mapYSup, N, xy);
486 // vectorPrint("y map", xy, nYc);
487 double yPrimitives[nYc];
488 axe = 1;
489 mathiesonPrimitive(xy, nYc, axe, chamberId, yPrimitives);
490
491 for (int i = 0; i < N; i++) {
492 Integrals[i] = (xPrimitives[mapXSup[i]] - xPrimitives[mapXInf[i]]) * (yPrimitives[mapYSup[i]] - yPrimitives[mapYInf[i]]);
493 // printf(" i=%d mapXInf=%d mapXSup=%d mapYInf=%d mapYSup=%d xyIntegrals=%f, %f \n", i,
494 // mapXInf[i], mapXSup[i], mapYInf[i], mapYSup[i], xPrimitives[mapXSup[i]] - xPrimitives[mapXInf[i]],
495 // yPrimitives[mapYSup[i]] - yPrimitives[mapYInf[i]]);
496 }
497
498 // vectorPrint("xPrimitives", xPrimitives, nXc);
499 // vectorPrint("yPrimitives", yPrimitives, nYc);
500
501 } else {
502
503 if (useSpline) {
504 double lBoundPrim[N], uBoundPrim[N], xIntegrals[N], yIntegrals[N];
505 int axe = 0;
506 // mathiesonPrimitive(xInf, N, axe, chamberId, lBoundPrim);
507 splineMathiesonPrimitive(xInf, N, axe, chamberId, lBoundPrim);
508 // mathiesonPrimitive(xSup, N, axe, chamberId, uBoundPrim);
509 splineMathiesonPrimitive(xSup, N, axe, chamberId, uBoundPrim);
510 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, xIntegrals);
511 // vectorPrint("xIntegrals analytics ", xIntegrals, N);
512 for (int i = 0; i < N; i++) {
513 if (xIntegrals[i] < 0.0) {
514 printf("??? %d x (%f %f) lInt=%f uInt%f xInt=%f\n", i, xInf[i], xSup[i], lBoundPrim[i], uBoundPrim[i], xIntegrals[i]);
515 throw std::out_of_range(
516 "[findLocalMaxWithPEM] ????");
517 }
518 }
519 axe = 1;
520 // mathiesonPrimitive(yInf, N, axe, chamberId, lBoundPrim);
521 splineMathiesonPrimitive(yInf, N, axe, chamberId, lBoundPrim);
522 // mathiesonPrimitive(ySup, N, axe, chamberId, uBoundPrim);
523 splineMathiesonPrimitive(ySup, N, axe, chamberId, uBoundPrim);
524 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, yIntegrals);
525 // vectorPrint("yIntegrals analytics ", yIntegrals, N);
526 vectorMultVector(xIntegrals, yIntegrals, N, Integrals);
527 // Invald ????
528
529 for (int i = 0; i < N; i++) {
530 if (yIntegrals[i] < 0.0) {
531 printf("??? %d y (%f %f) lInt=%f uInt%f yInt=%f\n", i, yInf[i], ySup[i], lBoundPrim[i], uBoundPrim[i], yIntegrals[i]);
532 throw std::out_of_range(
533 "[findLocalMaxWithPEM] ????");
534 }
535 } // vectorPrint("Integrals analytics", Integrals, N);
536
537 /* ??????????????????????
538 axe = 0;
539 splineMathiesonPrimitive( xInf, N, axe, chamberId, lBoundPrim );
540 // vectorPrint("x lBoundPrim spline ", lBoundPrim, N);
541 splineMathiesonPrimitive( xSup, N, axe, chamberId, uBoundPrim );
542 vectorAddVector( uBoundPrim, -1.0, lBoundPrim, N, xIntegrals);
543 // vectorPrint("xIntegrals spline", xIntegrals, N);
544 axe = 1;
545 splineMathiesonPrimitive( yInf, N, axe, chamberId, lBoundPrim );
546 splineMathiesonPrimitive( ySup, N, axe, chamberId, uBoundPrim );
547 vectorAddVector( uBoundPrim, -1.0, lBoundPrim, N, yIntegrals);
548
549 // vectorPrint("yIntegrals spline", yIntegrals, N);
550 */
551
552 vectorMultVector(xIntegrals, yIntegrals, N, Integrals);
553 // vectorPrint("Integrals spline", Integrals, N);
554
555 } else {
556 // Returning array: Charge Integral on all the pads
557 //
558 if (chamberId <= 2) {
559 mathiesonType = 0;
560 } else {
561 mathiesonType = 1;
562 }
563 //
564 // Select Mathieson coef.
565 double curK2x = K2x[mathiesonType];
566 double curK2y = K2y[mathiesonType];
567 double curSqrtK3x = sqrtK3x[mathiesonType];
568 double curSqrtK3y = sqrtK3y[mathiesonType];
569 double curK4x = K4x[mathiesonType];
570 double curK4y = K4y[mathiesonType];
571 double curInvPitch = invPitch[mathiesonType];
572 double cst2x = curK2x * curInvPitch;
573 double cst2y = curK2y * curInvPitch;
574 double cst4 = 4.0 * curK4x * curK4y;
575 double uInf, uSup, vInf, vSup;
576
577 for (int i = 0; i < N; i++) {
578 // x/u
579 uInf = curSqrtK3x * tanh(cst2x * xInf[i]);
580 uSup = curSqrtK3x * tanh(cst2x * xSup[i]);
581 // y/v
582 vInf = curSqrtK3y * tanh(cst2y * yInf[i]);
583 vSup = curSqrtK3y * tanh(cst2y * ySup[i]);
584 //
585 Integrals[i] = cst4 * (atan(uSup) - atan(uInf)) * (atan(vSup) - atan(vInf));
586 // printf(" Ix=%10.6g Iy=%10.6g\n", 2*curK4x * (atan(uSup) - atan(uInf)), 2*curK4y * (atan(vSup) - atan(vInf)));
587 // printf(" xyInfSup %2d [%10.6g, %10.6g] x [%10.6g, %10.6g]-> %10.6g * %10.6g = %10.6g\n",
588 // i, xInf[i], xSup[i], yInf[i], ySup[i], Integrals[i], 2.0 * curK4x*(atan(uSup) - atan(uInf)), 2.0 * curK4y*(atan(vSup) - atan(vInf)) ) ;
589 }
590 // printf(" I[0..%3ld] = %f, %f, ... %f\n", N-1, Integrals[0], Integrals[1],
591 // Integrals[N-1]);
592 }
593 }
594 // CHECK
596 checkIntegrals(xInf, xSup, yInf, ySup, Integrals, chamberId, N);
597 }
598}
599
600void compute2DMathiesonMixturePadIntegrals(const double* xyInfSup0,
601 const double* theta, int N, int K,
602 int chamberId, double Integrals[])
603{
604 // Returning array: Charge Integral on all the pads
605 // Remarks:
606 // - This fct is a cumulative one, as a result it should be set to zero
607 // before calling it
608 vectorSetZero(Integrals, N);
609 const double* xInf0 = getConstXInf(xyInfSup0, N);
610 const double* yInf0 = getConstYInf(xyInfSup0, N);
611 const double* xSup0 = getConstXSup(xyInfSup0, N);
612 const double* ySup0 = getConstYSup(xyInfSup0, N);
613 //
614 const double* muX = getConstMuX(theta, K);
615 const double* muY = getConstMuY(theta, K);
616 const double* w = getConstW(theta, K);
617
618 double z[N];
619 double xyInfSup[4 * N];
620 double* xInf = getXInf(xyInfSup, N);
621 double* yInf = getYInf(xyInfSup, N);
622 double* xSup = getXSup(xyInfSup, N);
623 double* ySup = getYSup(xyInfSup, N);
624 for (int k = 0; k < K; k++) {
625 vectorAddScalar(xInf0, -muX[k], N, xInf);
626 vectorAddScalar(xSup0, -muX[k], N, xSup);
627 vectorAddScalar(yInf0, -muY[k], N, yInf);
628 vectorAddScalar(ySup0, -muY[k], N, ySup);
629 compute2DPadIntegrals(xInf, xSup, yInf, ySup, N, chamberId, z);
630 // printf("Vector Sum %g\n", vectorSum(z, N) );
631 vectorAddVector(Integrals, w[k], z, N, Integrals);
632 }
633}
634
635bool checkIntegrals(const double* xInf, const double* xSup, const double* yInf, const double* ySup,
636 const double* integralsToCheck, int chId, int N)
637{
638 double lBoundPrim[N], uBoundPrim[N];
639 double xIntegrals[N], yIntegrals[N], Integrals[N];
640 // ??? find the reason for high value
641 double precision = 5.e-5;
642 int axe = 0;
643 mathiesonPrimitive(xInf, N, axe, chId, lBoundPrim);
644 mathiesonPrimitive(xSup, N, axe, chId, uBoundPrim);
645 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, xIntegrals);
646 /*
647 for (int i=0; i < N; i++) {
648 if ( xIntegrals[i] >= 0.0) {
649 printf("i=%d xInf=%f xSup=%f, uBoundPrim=%f lBoundPrim=%f\n", i,
650 xInf[i], xSup[i], uBoundPrim[i], lBoundPrim[i]);
651 }
652 }
653 */
654 axe = 1;
655 mathiesonPrimitive(yInf, N, axe, chId, lBoundPrim);
656 mathiesonPrimitive(ySup, N, axe, chId, uBoundPrim);
657 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, yIntegrals);
658 // vectorPrint("yIntegrals analytics ", yIntegrals, N);
659 vectorMultVector(xIntegrals, yIntegrals, N, Integrals);
660 bool ok = true;
661 for (int i = 0; i < N; i++) {
662 if (std::fabs(integralsToCheck[i] - Integrals[i]) > precision) {
663 printf("i=%d xInf=%f xSup=%f, yInf=%f ySup=%f, reference=%f check value=%f\n", i,
664 xInf[i], xSup[i], yInf[i], ySup[i], Integrals[i], integralsToCheck[i]);
665 ok = false;
666 throw std::out_of_range("[checkIntegral] bad integral value");
667 }
668 }
669
670 return ok;
671}
672
673void computeFastCij(const Pads& pads, const Pads& pixel, double Cij[])
674{
675 // Compute the Charge Integral Cij of pads (j index), considering the
676 // center of the Mathieson fct on a pixel (i index)
677 // Use the fact that the charge integral CI(x,y) = CI(x) * CI(y)
678 // to reduce the computation cost
679 // CI(x) is store in PadIntegralX
680 // CI(y) is store in PadIntegralY
681 // A sub-sampling of CI(x_i + k*minDx) (or CI(y_i + l*minDY)) is used
682 // by taking the mininimun of pads.dx(pads.dy) to discretize the x/y space
683 //
684 // CI(x)/CI(y) are computed if they are requested.
685 //
686 // Returning array: Charge Integral on all the pads Cij[]
687
689 printf(
690 "[computeFastCij] exception: bad representation (mode) of pads in "
691 "computeCij (padMode=%d, pixelMode=%d)\n",
692 (int)pads.mode, (int)pixel.mode);
693 throw std::overflow_error("Bad mode");
694 return;
695 }
696 int N = pads.getNbrOfPads();
697 int K = pixel.getNbrOfPads();
698 // Pads
699 int chId = pads.getChamberId();
700 const double* xInf0 = pads.getXInf();
701 const double* yInf0 = pads.getYInf();
702 const double* xSup0 = pads.getXSup();
703 const double* ySup0 = pads.getYSup();
704 // Pixels
705 const double* muX = pixel.getX();
706 const double* muY = pixel.getY();
707
708 double zInf[N];
709 double zSup[N];
710 int axe;
711
712 // Loop on Pixels
713 std::map<int, double*> xMap;
714 std::map<int, double*> yMap;
715 for (int k = 0; k < K; k++) {
716 // Calculate the indexes in the 1D charge integral
717 // Error on pad position > 10-3 cm
718 int xCode = (int)(muX[k] * 1000 + 0.5);
719 int yCode = (int)(muY[k] * 1000 + 0.5);
720 if (xMap.find(xCode) == xMap.end()) {
721 // Not yet computed
722 vectorAddScalar(xInf0, -muX[k], N, zInf);
723 vectorAddScalar(xSup0, -muX[k], N, zSup);
724 axe = 0;
725 double* xIntegrals = new double[N];
726 compute1DPadIntegrals(zInf, zSup, N, axe, chId, xIntegrals);
727 xMap[xCode] = xIntegrals;
728 }
729 if (yMap.find(yCode) == yMap.end()) {
730 // Not yet computed
731 vectorAddScalar(yInf0, -muY[k], N, zInf);
732 vectorAddScalar(ySup0, -muY[k], N, zSup);
733 axe = 1;
734 double* yIntegrals = new double[N];
735 compute1DPadIntegrals(zInf, zSup, N, axe, chId, yIntegrals);
736 yMap[yCode] = yIntegrals;
737 }
738 // Compute IC(xy) = IC(x) * IC(y)
739 vectorMultVector(xMap[xCode], yMap[yCode], N, &Cij[N * k]);
740 //
741 // Check
743 double xInf[N], xSup[N];
744 double yInf[N], ySup[N];
745 double lBoundPrim[N], uBoundPrim[N];
746 double xIntegrals[N], yIntegrals[N], Integrals[N];
747 // printf("pad xyPad[0]= %f %f \n", (xSup0[0] - xInf0[0])*0.5, (ySup0[0] - yInf0[0])*0.5);
748 // printf("pad xyPad[0]= %f %f \n", xSup0[0], ySup0[0]);
749 // printf("pad xyPix[0]= %f %f \n", muX[k], muY[k]);
750 vectorAddScalar(xInf0, -muX[k], N, xInf);
751 vectorAddScalar(xSup0, -muX[k], N, xSup);
752 vectorAddScalar(yInf0, -muY[k], N, yInf);
753 vectorAddScalar(ySup0, -muY[k], N, ySup);
754 checkIntegrals(xInf, xSup, yInf, ySup, &Cij[N * k], chId, N);
755 }
756 }
757 // Free map
758 for (auto it = xMap.begin(); it != xMap.end(); ++it) {
759 delete[] it->second;
760 }
761 for (auto it = yMap.begin(); it != yMap.end(); ++it) {
762 delete[] it->second;
763 }
764}
765
766void computeFastCijV0(const Pads& pads, const Pads& pixel, double Cij[])
767{
768 // Compute the Charge Integral Cij of pads (j index), considering the
769 // center of the Mathieson fct on a pixel (i index)
770 // Use the fact that the charge integral CI(x,y) = CI(x) * CI(y)
771 // to reduce the computation cost
772 // CI(x) is store in PadIntegralX
773 // CI(y) is store in PadIntegralY
774 // A sub-sampling of CI(x_i + k*minDx) (or CI(y_i + l*minDY)) is used
775 // by taking the mininimun of pads.dx(pads.dy) to discretize the x/y space
776 //
777 // CI(x)/CI(y) are computed if they are requested.
778 //
779 // Returning array: Charge Integral on all the pads Cij[]
780
782 printf(
783 "[computeFastCij] exception: bad representation (mode) of pads in "
784 "computeCij (padMode=%d, pixelMode=%d)\n",
785 (int)pads.mode, (int)pixel.mode);
786 throw std::overflow_error("Bad mode");
787 return;
788 }
789 int N = pads.getNbrOfPads();
790 int K = pixel.getNbrOfPads();
791 // Pads
792 int chId = pads.getChamberId();
793 const double* xInf0 = pads.getXInf();
794 const double* yInf0 = pads.getYInf();
795 const double* xSup0 = pads.getXSup();
796 const double* ySup0 = pads.getYSup();
797 // Pixels
798 const double* muX = pixel.getX();
799 const double* muY = pixel.getY();
800 //
801 double xPixMin = vectorMin(muX, K);
802 double xPixMax = vectorMax(muX, K);
803 double yPixMin = vectorMin(muY, K);
804 double yPixMax = vectorMax(muY, K);
805 double dxMinPix = vectorMin(pixel.getDX(), K);
806 double dyMinPix = vectorMin(pixel.getDY(), K);
807 // Sampling of PadIntegralX/PadIntegralY
808 int nXPixels = (int)((xPixMax - xPixMin) / dxMinPix + 0.5) + 1;
809 int nYPixels = (int)((yPixMax - yPixMin) / dyMinPix + 0.5) + 1;
810 //
811 // PadIntegralX/PadIntegralY allocation and init with -1
812 // ??? double PadIntegralX[nXPixels][N];
813 // double PadIntegralY[nYPixels][N];
814 double* PadIntegralX = new double[nXPixels * N];
815 double* PadIntegralY = new double[nYPixels * N];
816 // Inv. printf("??? nXPixels=%d, xPixMin=%f, xPixMax=%f, dxMinPix=%f, nPads=%d\n", nXPixels, xPixMin, xPixMax, dxMinPix, N);
817 // printf("??? nYPixels=%d, yPixMin=%f, yPixMax=%f, dyMinPix=%f, nPads=%d\n", nYPixels, yPixMin, yPixMax, dyMinPix, N);
818 vectorSet((double*)PadIntegralX, -1.0, nXPixels * N);
819 vectorSet((double*)PadIntegralY, -1.0, nYPixels * N);
820 double zInf[N];
821 double zSup[N];
822 int axe;
823 /*
824 for (int kx=0; kx < nXPixels; kx++) {
825 double x = xPixMin + kx * dxPix;
826 vectorAddScalar( xInf0, - x, N, zInf );
827 vectorAddScalar( xSup0, - x, N, zSup );
828 compute1DPadIntegrals( zInf, zSup, N, chId, xAxe, PadIntegralX[kx] );
829 }
830 xAxe = false;
831 for (int ky=0; ky < nYPixels; ky++) {
832 double y = yPixMin + ky * dyPix;
833 vectorAddScalar( yInf0, - y, N, zInf );
834 vectorAddScalar( ySup0, - y, N, zSup );
835 compute1DPadIntegrals( zInf, zSup, N, chId, xAxe, PadIntegralY[ky] );
836 }
837 */
838
839 // Loop on Pixels
840 for (int k = 0; k < K; k++) {
841 // Calculate the indexes in the 1D charge integral
842 // PadIntegralX:PadIntegralY
843 int xIdx = (int)((muX[k] - xPixMin) / dxMinPix + 0.5);
844 int yIdx = (int)((muY[k] - yPixMin) / dyMinPix + 0.5);
845 // compute2DPadIntegrals( xInf, xSup, yInf, ySup, N, chId, &Cij[N*k] );
846 // Cij[ N*k + p] = PadIntegralX( k, xIdx) * PadIntegralY( k, yIdx);
847 // printf("k=%d, mu[k]=(%f, %f) Sum_pads Ck = %g\n", k, muX[k], muY[k],
848 // vectorSum( &Cij[N*k], N) );
849 if (PadIntegralX[xIdx * N + 0] == -1) {
850 // Not yet computed
851 vectorAddScalar(xInf0, -muX[k], N, zInf);
852 vectorAddScalar(xSup0, -muX[k], N, zSup);
853 axe = 0;
854 compute1DPadIntegrals(zInf, zSup, N, axe, chId, &PadIntegralX[xIdx * N + 0]);
855 }
856 if (PadIntegralY[yIdx * N + 0] == -1) {
857 // Not yet computed
858 vectorAddScalar(yInf0, -muY[k], N, zInf);
859 vectorAddScalar(ySup0, -muY[k], N, zSup);
860 axe = 1;
861 compute1DPadIntegrals(zInf, zSup, N, axe, chId, &PadIntegralY[yIdx * N + 0]);
862 }
863 // Compute IC(xy) = IC(x) * IC(y)
864 vectorMultVector(&PadIntegralX[xIdx * N + 0], &PadIntegralY[yIdx * N + 0], N, &Cij[N * k]);
865
866 double xInf[N], xSup[N];
867 double yInf[N], ySup[N];
868 double lBoundPrim[N], uBoundPrim[N];
869 double xIntegrals[N], yIntegrals[N], Integrals[N];
870
871 vectorAddScalar(xInf0, -muX[k], N, xInf);
872 vectorAddScalar(xSup0, -muX[k], N, xSup);
873 vectorAddScalar(yInf0, -muY[k], N, yInf);
874 vectorAddScalar(ySup0, -muY[k], N, ySup);
875 double integral;
876 int axe = 0;
877 mathiesonPrimitive(xInf, N, axe, chId, lBoundPrim);
878 mathiesonPrimitive(xSup, N, axe, chId, uBoundPrim);
879 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, xIntegrals);
880 axe = 1;
881 mathiesonPrimitive(yInf, N, axe, chId, lBoundPrim);
882 mathiesonPrimitive(ySup, N, axe, chId, uBoundPrim);
883 vectorAddVector(uBoundPrim, -1.0, lBoundPrim, N, yIntegrals);
884 // vectorPrint("yIntegrals analytics ", yIntegrals, N);
885 vectorMultVector(xIntegrals, yIntegrals, N, Integrals);
886 for (int i = 0; i < N; i++) {
887 // compute2DPadIntegrals(xInf[i], xSup, yInf, ySup, 1, chId, &integral);
888 if (std::fabs(Cij[N * k + i] - Integrals[i]) > 1.0e-6) {
889 printf("i(pixel)=%d j(pad)=%d cij=%f xInt=%f yInt=%f fastcij=%f xFast=%f yFast=%f\n", k, i,
890 Integrals[i], xIntegrals[i], yIntegrals[i], Cij[N * k + i], PadIntegralX[xIdx * N + i], PadIntegralY[yIdx * N + i]);
891 }
892 }
893 }
894 delete[] PadIntegralX;
895 delete[] PadIntegralY;
896}
897
898void computeCij(const Pads& pads, const Pads& pixel, double Cij[])
899{
900 // Compute the Charge Integral Cij of pads (j index), considering the
901 // center of the Mathieson fct on a pixel (i index)
902 //
903 // Returning array: Charge Integral on all the pads Cij[]
904
906 printf(
907 "computeFastCij] exception: bad representation (mode) of pads in "
908 "computeCij (padMode=%d, pixelMode=%d)\n",
909 (int)pads.mode, (int)pixel.mode);
910 throw std::overflow_error("Bad mode");
911 return;
912 }
913 int N = pads.getNbrOfPads();
914 int K = pixel.getNbrOfPads();
915 int chId = pads.getChamberId();
916 const double* xInf0 = pads.getXInf();
917 const double* yInf0 = pads.getYInf();
918 const double* xSup0 = pads.getXSup();
919 const double* ySup0 = pads.getYSup();
920
921 //
922 const double* muX = pixel.getX();
923 const double* muY = pixel.getY();
924
925 double xInf[N];
926 double yInf[N];
927 double xSup[N];
928 double ySup[N];
929
930 for (int k = 0; k < K; k++) {
931 vectorAddScalar(xInf0, -muX[k], N, xInf);
932 vectorAddScalar(xSup0, -muX[k], N, xSup);
933 vectorAddScalar(yInf0, -muY[k], N, yInf);
934 vectorAddScalar(ySup0, -muY[k], N, ySup);
935 compute2DPadIntegrals(xInf, xSup, yInf, ySup, N, chId, &Cij[N * k]);
936 // printf("k=%d, mu[k]=(%f, %f) Sum_pads Ck = %g\n", k, muX[k], muY[k],
937 // vectorSum( &Cij[N*k], N) );
938 }
939}
940
941void checkCij(const Pads& pads, const Pads& pixels, const double* checkCij, int mode)
942{
943 // Mode : 0 (nothing), 1 (info), 2 (detail), -1 (exception)
944 int nPads = pads.getNbrOfPads();
945 int nPixels = pixels.getNbrOfPads();
946 double* Cij = new double[nPads * nPixels];
947 double* diffCij = new double[nPads * nPixels];
948 double precision = 2.0e-5;
949 computeCij(pads, pixels, Cij);
950 vectorAddVector(Cij, -1, checkCij, nPads * nPixels, diffCij);
951 vectorAbs(diffCij, nPads * nPixels, diffCij);
952 double minDiff = vectorMin(diffCij, nPads * nPixels);
953 double maxDiff = vectorMax(diffCij, nPads * nPixels);
954 int argMax = vectorArgMax(diffCij, nPads * nPixels);
955 // printf("\n\n nPads, nPixels %d %d\n", nPads, nPixels);
956 int iIdx = argMax / nPads;
957 int jIdx = argMax % nPads;
958 if ((maxDiff > precision) && (mode != 0)) {
959 printf("\n\n[checkCij] min/max(checkCij-Cij)=(%f, %f) argmin/max=(i=%d, j=%d)\n",
960 minDiff, maxDiff, iIdx, jIdx);
961 printf("\n checkCij=%f differ from Cij=%f\n", checkCij[iIdx * nPads + jIdx], Cij[iIdx * nPads + jIdx]);
962 }
963
964 if ((maxDiff > precision) && (mode > 1)) {
965 for (int k = 0; k < nPixels; k++) {
966 for (int l = 0; l < nPads; l++) {
967 if (diffCij[k * nPads + l] > precision) {
968 printf("pad=%d pixel=%d checkCij=%f Cij=%f diff=%f\n", l, k, checkCij[k * nPads + l], Cij[k * nPads + l], diffCij[k * nPads + l]);
969 }
970 }
971 }
972 // printf("findLocalMaxWithPEM: WARNING maxDiff(Cij)=%f\n", maxDiff);
973 }
974 if ((maxDiff > precision) && (mode == -1)) {
975 throw std::out_of_range("[checkCij] bad Cij value");
976 }
977 delete[] Cij;
978 delete[] diffCij;
979}
980
981// theta
982double* getVarX(double* theta, int K) { return &theta[0 * K]; };
983double* getVarY(double* theta, int K) { return &theta[1 * K]; };
984double* getMuX(double* theta, int K) { return &theta[2 * K]; };
985double* getMuY(double* theta, int K) { return &theta[3 * K]; };
986double* getW(double* theta, int K) { return &theta[4 * K]; };
987double* getMuAndW(double* theta, int K) { return &theta[2 * K]; };
988//
989const double* getConstVarX(const double* theta, int K)
990{
991 return &theta[0 * K];
992};
993const double* getConstVarY(const double* theta, int K)
994{
995 return &theta[1 * K];
996};
997const double* getConstMuX(const double* theta, int K) { return &theta[2 * K]; };
998const double* getConstMuY(const double* theta, int K) { return &theta[3 * K]; };
999const double* getConstW(const double* theta, int K) { return &theta[4 * K]; };
1000const double* getConstMuAndW(const double* theta, int K)
1001{
1002 return &theta[2 * K];
1003};
1004
1005// xyDxy
1006double* getX(double* xyDxy, int N) { return &xyDxy[0 * N]; };
1007double* getY(double* xyDxy, int N) { return &xyDxy[1 * N]; };
1008double* getDX(double* xyDxy, int N) { return &xyDxy[2 * N]; };
1009double* getDY(double* xyDxy, int N) { return &xyDxy[3 * N]; };
1010//
1011const double* getConstX(const double* xyDxy, int N) { return &xyDxy[0 * N]; };
1012const double* getConstY(const double* xyDxy, int N) { return &xyDxy[1 * N]; };
1013const double* getConstDX(const double* xyDxy, int N) { return &xyDxy[2 * N]; };
1014const double* getConstDY(const double* xyDxy, int N) { return &xyDxy[3 * N]; };
1015
1016// xySupInf
1017double* getXInf(double* xyInfSup, int N) { return &xyInfSup[0 * N]; };
1018double* getYInf(double* xyInfSup, int N) { return &xyInfSup[1 * N]; };
1019double* getXSup(double* xyInfSup, int N) { return &xyInfSup[2 * N]; };
1020double* getYSup(double* xyInfSup, int N) { return &xyInfSup[3 * N]; };
1021const double* getConstXInf(const double* xyInfSup, int N)
1022{
1023 return &xyInfSup[0 * N];
1024};
1025const double* getConstYInf(const double* xyInfSup, int N)
1026{
1027 return &xyInfSup[1 * N];
1028};
1029const double* getConstXSup(const double* xyInfSup, int N)
1030{
1031 return &xyInfSup[2 * N];
1032};
1033const double* getConstYSup(const double* xyInfSup, int N)
1034{
1035 return &xyInfSup[3 * N];
1036};
1037
1038void copyTheta(const double* theta0, int K0, double* theta, int K1, int K)
1039{
1040 const double* w = getConstW(theta0, K0);
1041 const double* muX = getConstMuX(theta0, K0);
1042 const double* muY = getConstMuY(theta0, K0);
1043 const double* varX = getConstVarX(theta0, K0);
1044 const double* varY = getConstVarY(theta0, K0);
1045 double* wm = getW(theta, K1);
1046 double* muXm = getMuX(theta, K1);
1047 double* muYm = getMuY(theta, K1);
1048 double* varXm = getVarX(theta, K1);
1049 double* varYm = getVarY(theta, K1);
1050 vectorCopy(w, K, wm);
1051 vectorCopy(muX, K, muXm);
1052 vectorCopy(muY, K, muYm);
1053 vectorCopy(varX, K, varXm);
1054 vectorCopy(varY, K, varYm);
1055}
1056
1057void copyXYdXY(const double* xyDxy0, int N0, double* xyDxy, int N1, int N)
1058{
1059 const double* X0 = getConstX(xyDxy0, N0);
1060 const double* Y0 = getConstY(xyDxy0, N0);
1061 const double* DX0 = getConstDX(xyDxy0, N0);
1062 const double* DY0 = getConstDY(xyDxy0, N0);
1063
1064 double* X = getX(xyDxy, N1);
1065 double* Y = getY(xyDxy, N1);
1066 double* DX = getDX(xyDxy, N1);
1067 double* DY = getDY(xyDxy, N1);
1068
1069 vectorCopy(X0, N, X);
1070 vectorCopy(Y0, N, Y);
1071 vectorCopy(DX0, N, DX);
1072 vectorCopy(DY0, N, DY);
1073}
1074
1075void printTheta(const char* str, double meanCharge, const double* theta, int K)
1076{
1077 const double* varX = getConstVarX(theta, K);
1078 const double* varY = getConstVarY(theta, K);
1079 const double* muX = getConstMuX(theta, K);
1080 const double* muY = getConstMuY(theta, K);
1081 const double* w = getConstW(theta, K);
1082
1083 printf("%s \n", str);
1084 printf(" k charge w muX muY sigX sigY\n");
1085 for (int k = 0; k < K; k++) {
1086 printf(" %.2d-th: %8.2g %6.3g %8.3g %8.3g %8.3g %8.3g\n", k, w[k] * meanCharge, w[k], muX[k],
1087 muY[k], sqrt(varX[k]), sqrt(varY[k]));
1088 }
1089}
1090void xyDxyToxyInfSup(const double* xyDxy, int nxyDxy, double* xyInfSup)
1091{
1092 const double* X = getConstX(xyDxy, nxyDxy);
1093 const double* Y = getConstY(xyDxy, nxyDxy);
1094 const double* DX = getConstDX(xyDxy, nxyDxy);
1095 const double* DY = getConstDY(xyDxy, nxyDxy);
1096 double* XInf = getXInf(xyInfSup, nxyDxy);
1097 double* XSup = getXSup(xyInfSup, nxyDxy);
1098 double* YInf = getYInf(xyInfSup, nxyDxy);
1099 double* YSup = getYSup(xyInfSup, nxyDxy);
1100 for (int k = 0; k < nxyDxy; k++) {
1101 // To avoid overwritting
1102 double xInf = X[k] - DX[k];
1103 double xSup = X[k] + DX[k];
1104 double yInf = Y[k] - DY[k];
1105 double ySup = Y[k] + DY[k];
1106 //
1107 XInf[k] = xInf;
1108 YInf[k] = yInf;
1109 XSup[k] = xSup;
1110 YSup[k] = ySup;
1111 }
1112}
1113
1114// Mask operations
1115void maskedCopyXYdXY(const double* xyDxy, int nxyDxy, const Mask_t* mask,
1116 int nMask, double* xyDxyMasked, int nxyDxyMasked)
1117{
1118 const double* X = getConstX(xyDxy, nxyDxy);
1119 const double* Y = getConstY(xyDxy, nxyDxy);
1120 const double* DX = getConstDX(xyDxy, nxyDxy);
1121 const double* DY = getConstDY(xyDxy, nxyDxy);
1122 double* Xm = getX(xyDxyMasked, nxyDxyMasked);
1123 double* Ym = getY(xyDxyMasked, nxyDxyMasked);
1124 double* DXm = getDX(xyDxyMasked, nxyDxyMasked);
1125 double* DYm = getDY(xyDxyMasked, nxyDxyMasked);
1126 vectorGather(X, mask, nMask, Xm);
1127 vectorGather(Y, mask, nMask, Ym);
1128 vectorGather(DX, mask, nMask, DXm);
1129 vectorGather(DY, mask, nMask, DYm);
1130}
1131
1132void maskedCopyToXYInfSup(const double* xyDxy, int ndxyDxy, const Mask_t* mask,
1133 int nMask, double* xyDxyMasked, int ndxyDxyMasked)
1134{
1135 const double* X = getConstX(xyDxy, ndxyDxy);
1136 const double* Y = getConstY(xyDxy, ndxyDxy);
1137 const double* DX = getConstDX(xyDxy, ndxyDxy);
1138 const double* DY = getConstDY(xyDxy, ndxyDxy);
1139 double* Xm = getX(xyDxyMasked, ndxyDxyMasked);
1140 double* Ym = getY(xyDxyMasked, ndxyDxyMasked);
1141 double* DXm = getDX(xyDxyMasked, ndxyDxyMasked);
1142 double* DYm = getDY(xyDxyMasked, ndxyDxyMasked);
1143 double* XmInf = getXInf(xyDxyMasked, ndxyDxyMasked);
1144 double* XmSup = getXSup(xyDxyMasked, ndxyDxyMasked);
1145 double* YmInf = getYInf(xyDxyMasked, ndxyDxyMasked);
1146 double* YmSup = getYSup(xyDxyMasked, ndxyDxyMasked);
1147 vectorGather(X, mask, nMask, Xm);
1148 vectorGather(Y, mask, nMask, Ym);
1149 vectorGather(DX, mask, nMask, DXm);
1150 vectorGather(DY, mask, nMask, DYm);
1151 for (int k = 0; k < nMask; k++) {
1152 // To avoid overwritting
1153 double xInf = Xm[k] - DXm[k];
1154 double xSup = Xm[k] + DXm[k];
1155 double yInf = Ym[k] - DYm[k];
1156 double ySup = Ym[k] + DYm[k];
1157 //
1158 XmInf[k] = xInf;
1159 YmInf[k] = yInf;
1160 XmSup[k] = xSup;
1161 YmSup[k] = ySup;
1162 }
1163}
1164
1165void maskedCopyTheta(const double* theta, int K, const Mask_t* mask, int nMask,
1166 double* maskedTheta, int maskedK)
1167{
1168 const double* w = getConstW(theta, K);
1169 const double* muX = getConstMuX(theta, K);
1170 const double* muY = getConstMuY(theta, K);
1171 const double* varX = getConstVarX(theta, K);
1172 const double* varY = getConstVarY(theta, K);
1173 double* wm = getW(maskedTheta, maskedK);
1174 double* muXm = getMuX(maskedTheta, maskedK);
1175 double* muYm = getMuY(maskedTheta, maskedK);
1176 double* varXm = getVarX(maskedTheta, maskedK);
1177 double* varYm = getVarY(maskedTheta, maskedK);
1178 vectorGather(w, mask, nMask, wm);
1179 vectorGather(muX, mask, nMask, muXm);
1180 vectorGather(muY, mask, nMask, muYm);
1181 vectorGather(varX, mask, nMask, varXm);
1182 vectorGather(varY, mask, nMask, varYm);
1183}
1184
1185void printXYdXY(const char* str, const double* xyDxy, int NMax, int N,
1186 const double* val1, const double* val2)
1187{
1188 const double* X = getConstX(xyDxy, NMax);
1189 const double* Y = getConstY(xyDxy, NMax);
1190 const double* DX = getConstDX(xyDxy, NMax);
1191 const double* DY = getConstDY(xyDxy, NMax);
1192
1193 printf("%s\n", str);
1194 int nPrint = 0;
1195 if (val1 != nullptr) {
1196 nPrint++;
1197 }
1198 if (val2 != nullptr) {
1199 nPrint++;
1200 }
1201 if ((nPrint == 1) && (val2 != nullptr)) {
1202 val1 = val2;
1203 }
1204
1205 if (nPrint == 0) {
1206 for (PadIdx_t i = 0; i < N; i++) {
1207 printf(" pad %2d: %9.3g %9.3g %9.3g %9.3g \n", i, X[i], Y[i], DX[i],
1208 DY[i]);
1209 }
1210 } else if (nPrint == 1) {
1211 for (PadIdx_t i = 0; i < N; i++) {
1212 printf(" pad %2d: %9.3g %9.3g %9.3g %9.3g - %9.3g \n", i, X[i], Y[i],
1213 DX[i], DY[i], val1[i]);
1214 }
1215 } else {
1216 for (PadIdx_t i = 0; i < N; i++) {
1217 printf(" pad %d: %9.3g %9.3g %9.3g %9.3g - %9.3g %9.3g \n", i, X[i],
1218 Y[i], DX[i], DY[i], val1[i], val2[i]);
1219 }
1220 }
1221}
1222
1223} // namespace mch
1224} // namespace o2
1225
1226// C Wrapper
1231
1232void o2_mch_compute2DPadIntegrals(const double* xInf, const double* xSup,
1233 const double* yInf, const double* ySup, int N,
1234 int chamberId, double Integrals[])
1235{
1236 o2::mch::compute2DPadIntegrals(xInf, xSup, yInf, ySup, N, chamberId,
1237 Integrals);
1238}
1239
1240void o2_mch_computeCij(const double* xyInfSup0, const double* pixel, int N,
1241 int K, int chamberId, double Cij[])
1242{
1243 // Returning array: Charge Integral on all the pads
1244 // Remarks:
1245 // - This fct is a cumulative one, as a result it should be set to zero
1246 // before calling it
1247 const double* xInf0 = o2::mch::getConstXInf(xyInfSup0, N);
1248 const double* yInf0 = o2::mch::getConstYInf(xyInfSup0, N);
1249 const double* xSup0 = o2::mch::getConstXSup(xyInfSup0, N);
1250 const double* ySup0 = o2::mch::getConstYSup(xyInfSup0, N);
1251 //
1252 const double* muX = o2::mch::getConstMuX(pixel, K);
1253 const double* muY = o2::mch::getConstMuY(pixel, K);
1254 const double* w = o2::mch::getConstW(pixel, K);
1255
1256 double z[N];
1257 double xyInfSup[4 * N];
1258 double* xInf = o2::mch::getXInf(xyInfSup, N);
1259 double* yInf = o2::mch::getYInf(xyInfSup, N);
1260 double* xSup = o2::mch::getXSup(xyInfSup, N);
1261 double* ySup = o2::mch::getYSup(xyInfSup, N);
1262 for (int k = 0; k < K; k++) {
1263 o2::mch::vectorAddScalar(xInf0, -muX[k], N, xInf);
1264 o2::mch::vectorAddScalar(xSup0, -muX[k], N, xSup);
1265 o2::mch::vectorAddScalar(yInf0, -muY[k], N, yInf);
1266 o2::mch::vectorAddScalar(ySup0, -muY[k], N, ySup);
1267 o2_mch_compute2DPadIntegrals(xInf, xSup, yInf, ySup, N, chamberId,
1268 &Cij[N * k]);
1269 // printf("Vector Sum %g\n", vectorSum(z, N) );
1270 }
1271}
1272
1274 const double* theta, int N,
1275 int K, int chamberId,
1276 double Integrals[])
1277{
1279 chamberId, Integrals);
1280}
Clustering and fifting parameters.
int32_t i
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
Class for time synchronization of RawReader instances.
const double * getX() const
Definition PadsPEM.h:91
PadMode mode
Definition PadsPEM.h:58
const double * getXSup() const
Definition PadsPEM.h:97
int getChamberId() const
Definition PadsPEM.h:106
@ xydxdyMode
x, y, dx, dy pad coordinates
@ xyInfSupMode
xInf=x, xSup=dx, yInf=y, ySup=dy pad coordinates
static int getNbrOfPads(const Pads *pads)
Definition PadsPEM.h:63
const double * getDY() const
Definition PadsPEM.h:94
const double * getYInf() const
Definition PadsPEM.h:96
const double * getY() const
Definition PadsPEM.h:92
const double * getXInf() const
Definition PadsPEM.h:95
const double * getYSup() const
Definition PadsPEM.h:98
const double * getDX() const
Definition PadsPEM.h:93
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum mode
Definition glcorearb.h:266
GLuint GLfloat GLfloat GLfloat x1
Definition glcorearb.h:5034
GLdouble f
Definition glcorearb.h:310
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLint GLint GLsizei GLint GLenum GLenum const void * pixels
Definition glcorearb.h:275
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLint GLuint mask
Definition glcorearb.h:291
GLenum GLint GLint * precision
Definition glcorearb.h:1899
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
void o2_mch_compute2DMathiesonMixturePadIntegrals(const double *xyInfSup0, const double *theta, int N, int K, int chamberId, double Integrals[])
void o2_mch_computeCij(const double *xyInfSup0, const double *pixel, int N, int K, int chamberId, double Cij[])
void o2_mch_initMathieson()
void o2_mch_compute2DPadIntegrals(const double *xInf, const double *xSup, const double *yInf, const double *ySup, int N, int chamberId, double Integrals[])
void initMathieson(int useSpline_, int useCache_)
Definition mathieson.cxx:73
const double sqrtK3x1_2
Definition mathieson.cxx:29
const double * getConstXInf(const double *xyInfSup, int N)
void maskedCopyToXYInfSup(const double *xyDxy, int ndxyDxy, const Mask_t *mask, int nMask, double *xyDxyMasked, int ndxyDxyMasked)
void maskedCopyTheta(const double *theta, int K, const Mask_t *mask, int nMask, double *maskedTheta, int maskedK)
double * getMuY(double *theta, int K)
const double * getConstMuAndW(const double *theta, int K)
void printTheta(const char *str, double meanCharge, const double *theta, int K)
const double sqrtK3y1_2
Definition mathieson.cxx:30
void compute2DPadIntegrals(const double *xInf, const double *xSup, const double *yInf, const double *ySup, int N, int chamberId, double Integrals[])
void printXYdXY(const char *str, const double *xyDxy, int NMax, int N, const double *val1, const double *val2)
double * getDX(double *xyDxy, int N)
bool checkIntegrals(const double *xInf, const double *xSup, const double *yInf, const double *ySup, const double *integralsToCheck, int chId, int N)
CompressedPads_t * compressPads(const double *xInf, const double *xSup, const double *yInf, const double *ySup, int N)
void computeFastCij(const Pads &pads, const Pads &pixel, double Cij[])
const double * getConstDX(const double *xyDxy, int N)
void splineMathiesonPrimitive(const double *x, int N, int axe, int chamberId, double *mPrimitive)
double * getX(double *xyDxy, int N)
double * getY(double *xyDxy, int N)
void compute1DPadIntegrals(const double *xyInf, const double *xySup, int N, double xy0, int axe, int chamberId, double *integrals)
double * getYSup(double *xyInfSup, int N)
void maskedCopyXYdXY(const double *xyDxy, int nxyDxy, const Mask_t *mask, int nMask, double *xyDxyMasked, int nxyDxyMasked)
const double sqrtK3y3_10
Definition mathieson.cxx:34
double * getYInf(double *xyInfSup, int N)
const double * getConstVarY(const double *theta, int K)
void computeCij(const Pads &pads, const Pads &pixel, double Cij[])
void copyTheta(const double *theta0, int K0, double *theta, int K1, int K)
double * getMuAndW(double *theta, int K)
int useSpline
Definition mathieson.cxx:48
const double sqrtK3x3_10
Definition mathieson.cxx:33
void initSplineMathiesonPrimitive()
Definition mathieson.cxx:94
SplineCoef * splineCoef[2][2]
Definition mathieson.cxx:49
const double * getConstW(const double *theta, int K)
double * getVarX(double *theta, int K)
void compute2DMathiesonMixturePadIntegrals(const double *xyInfSup0, const double *theta, int N, int K, int chamberId, double Integrals[])
int useCache
Definition mathieson.cxx:56
void computeCompressed2DPadIntegrals(CompressedPads_t *compressedPads, double xShift, double yShift, int N, int chamberId, double Integrals[])
double * splineXY
Definition mathieson.cxx:53
void checkCij(const Pads &pads, const Pads &pixels, const double *checkCij, int mode)
void xyDxyToxyInfSup(const double *xyDxy, int nxyDxy, double *xyInfSup)
const double * getConstDY(const double *xyDxy, int N)
double * getMuX(double *theta, int K)
void copyXYdXY(const double *xyDxy0, int N0, double *xyDxy, int N1, int N)
double * getXSup(double *xyInfSup, int N)
double * getW(double *theta, int K)
const double * getConstYInf(const double *xyInfSup, int N)
const double * getConstY(const double *xyDxy, int N)
const double pitch1_2
Definition mathieson.cxx:31
const double * getConstX(const double *xyDxy, int N)
const double * getConstMuX(const double *theta, int K)
void deleteCompressedPads(CompressedPads_t *compressedPads)
double * getVarY(double *theta, int K)
double * getDY(double *xyDxy, int N)
const double * getConstXSup(const double *xyInfSup, int N)
const double * getConstMuY(const double *theta, int K)
ClusterConfig clusterConfig
const double pitch3_10
Definition mathieson.cxx:35
int compressSameValues(const double *x1, const double *x2, int *map1, int *map2, int N, double *xCompress)
void compute1DMathieson(const double *xy, int N, int axe, int chamberId, double mathieson[])
const double * getConstVarX(const double *theta, int K)
int mathiesonType
Definition mathieson.cxx:37
void computeFastCijV0(const Pads &pads, const Pads &pixel, double Cij[])
void mathiesonPrimitive(const double *xy, int N, int axe, int chamberId, double mPrimitive[])
void computeSplineCoef(const double *xy, double xyStep, const double *f, int N, double leftDerivative, double rightDerivative, SplineCoef *splineCoef)
double * getXInf(double *xyInfSup, int N)
const double * getConstYSup(const double *xyInfSup, int N)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
const std::string str