SND@LHC Software
Loading...
Searching...
No Matches
ScifiCluster Namespace Reference

Classes

class  ScifiCluster
 Class of SciFi Tracker ######################################################. More...
 

Functions

 cm_to_channel (locpos, sipm_map=sipm_map, gaps_map=gaps_map, pitch=pitch, charr=charr, reverse=False, ch_max_num=ch_max_num, n_solid_sipms=n_solid_sipms)
 FUNCTIONS FOR SCIFICLASS ##################################################.
 
 channel_to_cm (channelpos, sipm_map=sipm_map, reverse=False, pitch=pitch)
 
 GetMatType (DetID)
 
 GetMatNum (DetID)
 
 GetMatLength ()
 
 GetMatQty ()
 
 GetSiPMs ()
 
 GetStationNum (DetID)
 
 global_to_local (DetID, globalpos)
 These functions assume that the detector is aligned with the beam ######.
 
 local_to_global (DetID, localpos)
 
 ly_loss_mean (distance, params=ly_loss_params)
 
 ly_attenuation (distance)
 
 cluster_width_mean (distance, params=cluster_width_mean_params)
 
 cluster_width_sigma (distance, params=cluster_width_sigma_params)
 
 cluster_width_random (distance, ly, percent=random_width_percent, cluster_width_min=cluster_width_min, cluster_width_max=cluster_width_max)
 
 approx_function (var, approx_data)
 
 edep_to_ly (energy, CDF_integral=CDF_integral, energy_range=energy_range, ly_linear_params=ly_linear_params, k_cdfs_corr=k_cdfs_corr, sigma_in_percent=sigma_in_percent, ly_CDF_params=ly_CDF_params, CDF_ly_params=CDF_ly_params, ly_CDF_landau_params=ly_CDF_landau_params)
 
 cluster_generator (amplitude, width, wmp, cluster_width_max=cluster_width_max, chpos_min=chpos_min, chpos_max=chpos_max)
 
 is_realistic (cluster, width)
 
 create_cluster (amplitude, width, wmp)
 
 fcheck_wall (cluster, wmp)
 
 weigthed_mean_pos (cluster)
 

Variables

int sipm_hor_pos = +1
 GLOBAL VALUES AND PARAMETERS ############.
 
int sipm_vert_pos = -1
 
bool is_sipm_hor_reversed = True
 
bool is_sipm_vert_reversed = True
 
float hw = 19.528
 
float pitch = 0.025
 
float charr = 1.6
 
float edge = 0.016
 
float gap = 0.022
 
float array = gap + 2.*charr
 
float biggap = 0.038
 
int ch_max_num = 128*12 - 1
 
int n_solid_sipms = 12
 
tuple sipm_map
 
dict gaps_map
 
float firstmatgap_startch = 521.5
 
float firstmatgap_endch = 523.5
 
float secmatgap_startch = 1046.5
 
float secmatgap_endch = 1048.5
 
float chpos_min = -0.5
 
float chpos_max = 1535.5
 
int cluster_width_min = 1
 
int cluster_width_max = 10
 
float ly_loss_params = 20.78, -0.26, 7.89, -3.06
 
float ly_loss_sigma_params = 6.8482, -0.5757, 3.4458
 
float cluster_width_mean_params = 0.6661, -0.006955, 2.163
 
float cluster_width_sigma_params = 1.103, -0.0005751
 
float random_width_percent = 0.01
 
float ly_min = 4.5
 
int ly_max = 104
 
float energy_range = 0.18, 0.477
 
float CDF_integral = 0.0185640424
 
float ly_linear_params = 332.882, -58.7085
 
float k_cdfs_corr = 0.993076766938
 
float sigma_in_percent = 0.01
 
int sigma_from_width = 1 / 4.
 
dict ly_CDF_params
 
dict ly_CDF_landau_params
 
dict CDF_ly_params
 
 snd_geo = ConfigRegistry.loadpy("$SNDSW_ROOT/geometry/shipLHC_geom_config.py")
 LOAD THE SCIFI GEOMETRY ############.
 
 scifimat_width = snd_geo.Scifi.scifimat_width
 
 scifimat_length = snd_geo.Scifi.scifimat_length
 
 scifimat_gap = snd_geo.Scifi.scifimat_gap
 
 scifi_x = snd_geo.Scifi.xdim
 
 scifi_y = snd_geo.Scifi.ydim
 
 scifimat_z = snd_geo.Scifi.scifimat_z
 
 nmats = snd_geo.Scifi.nmats
 
 n_planes = snd_geo.Scifi.nscifi
 
 n_sipms = snd_geo.Scifi.nsipms
 
float xshift = -27.8
 Need to automatize this part so it works for every geometry #########.
 
float yshift = +35.3
 

Function Documentation

◆ approx_function()

ScifiCluster.approx_function (   var,
  approx_data 
)
    This universal function substitutes the parameters to the function. 
    The parameters and the function are in the dictionary

Definition at line 466 of file ScifiCluster.py.

466def approx_function(var, approx_data):
467
468 """
469 This universal function substitutes the parameters to the function.
470 The parameters and the function are in the dictionary
471 """
472
473 for (left, right), (params, func) in approx_data.items():
474 if left <= var <= right:
475 return func(var, *params)
476 return False
477

◆ channel_to_cm()

ScifiCluster.channel_to_cm (   channelpos,
  sipm_map = sipm_map,
  reverse = False,
  pitch = pitch 
)
    It converts a particle position measured channels to a position measured 
    in cm. The SiPM map is used. The position is in the scifi module frame.

Definition at line 304 of file ScifiCluster.py.

304def channel_to_cm(channelpos, sipm_map=sipm_map, reverse=False, pitch=pitch):
305
306 """
307 It converts a particle position measured channels to a position measured
308 in cm. The SiPM map is used. The position is in the scifi module frame.
309 """
310
311 if reverse is True:
312 for left, _, ch_range in sipm_map:
313 if not ch_range is False:
314 ch_start, ch_end = ch_range
315 if ch_start <= channelpos <= ch_end:
316 return -(left + (channelpos - ch_start) * pitch)
317 if reverse is False:
318 for left, _, ch_range in sipm_map:
319 if not ch_range is False:
320 ch_start, ch_end = ch_range
321 if ch_start <= channelpos <= ch_end:
322 return left + (channelpos - ch_start) * pitch
323

◆ cluster_generator()

ScifiCluster.cluster_generator (   amplitude,
  width,
  wmp,
  cluster_width_max = cluster_width_max,
  chpos_min = chpos_min,
  chpos_max = chpos_max 
)
    It generates an event cluster with given weighted mean position in channels, width and 
    amplitude. 

    If right side of the cluster can be out of range, the maximum of the right side will be 
    right channel.

    At first an array [0, 0, 0, ... ] is generated which corresponds to the channels. 
    Next the cluster generated in the array. 
    Final array will be like [0, 0, ..., 1, 2, 5, 1, 0, ...], 
    [0, 17, 0, ...] or etc.

Definition at line 513 of file ScifiCluster.py.

514 chpos_min=chpos_min, chpos_max=chpos_max):
515
516 """
517 It generates an event cluster with given weighted mean position in channels, width and
518 amplitude.
519
520 If right side of the cluster can be out of range, the maximum of the right side will be
521 right channel.
522
523 At first an array [0, 0, 0, ... ] is generated which corresponds to the channels.
524 Next the cluster generated in the array.
525 Final array will be like [0, 0, ..., 1, 2, 5, 1, 0, ...],
526 [0, 17, 0, ...] or etc.
527 """
528
529 if int(ceil(wmp + 0.5 + cluster_width_max)) < chpos_max:
530 max_fired_ch = int(ceil(wmp + 0.5 + cluster_width_max))
531 else:
532 max_fired_ch = int(ceil(chpos_max))
533 cluster = [0 for _ in range(max_fired_ch)]
534 mean = wmp
535 if width == 1:
536 if wmp == chpos_min:
537 fired_channel = 0
538 elif wmp == chpos_max:
539 fired_channel = int(chpos_max)
540 else:
541 fired_channel = int(wmp + 0.5)
542 cluster[fired_channel] += amplitude
543 else:
544 sigma = width * sigma_from_width
545 for _ in range(amplitude):
546 fired_channel = int(round(random.gauss(mean, sigma)))
547 while not 0 <= fired_channel < len(cluster):
548 fired_channel = int(round(random.gauss(mean, sigma)))
549 cluster[fired_channel] += 1
550
551 return cluster
552

◆ cluster_width_mean()

ScifiCluster.cluster_width_mean (   distance,
  params = cluster_width_mean_params 
)
    It return a mean cluster width depending on the distance to SiPMs

Definition at line 428 of file ScifiCluster.py.

428def cluster_width_mean(distance, params=cluster_width_mean_params):
429
430 """
431 It return a mean cluster width depending on the distance to SiPMs
432 """
433
434 A, B, C = params
435 return exp(A + B*distance) + C
436

◆ cluster_width_random()

ScifiCluster.cluster_width_random (   distance,
  ly,
  percent = random_width_percent,
  cluster_width_min = cluster_width_min,
  cluster_width_max = cluster_width_max 
)
    It generates a cluster. The cluster have 'ly' photoelectrons. 
    The cluster width depends on the distance to SiPM

Definition at line 446 of file ScifiCluster.py.

447 cluster_width_min=cluster_width_min, cluster_width_max=cluster_width_max):
448
449 """
450 It generates a cluster. The cluster have 'ly' photoelectrons.
451 The cluster width depends on the distance to SiPM
452 """
453
454 mean = cluster_width_mean(distance)
455 sigma = cluster_width_sigma(distance)
456 if random.random() <= percent:
457 return random.randint(cluster_width_min, cluster_width_max)
458 random_width = int(round(random.gauss(mean - 1, sigma))) + 1
459
460 # Generate again if the width < minimal width and the light yield < the width
461 while random_width < 1 or ly < random_width:
462 random_width = int(round(random.gauss(mean - 1, sigma))) + 1
463
464 return random_width
465

◆ cluster_width_sigma()

ScifiCluster.cluster_width_sigma (   distance,
  params = cluster_width_sigma_params 
)
    It return a standard deviation of the mean cluster width depending on the distance to SiPMs

Definition at line 437 of file ScifiCluster.py.

437def cluster_width_sigma(distance, params=cluster_width_sigma_params):
438
439 """
440 It return a standard deviation of the mean cluster width depending on the distance to SiPMs
441 """
442
443 A, B = params
444 return A + B*distance
445

◆ cm_to_channel()

ScifiCluster.cm_to_channel (   locpos,
  sipm_map = sipm_map,
  gaps_map = gaps_map,
  pitch = pitch,
  charr = charr,
  reverse = False,
  ch_max_num = ch_max_num,
  n_solid_sipms = n_solid_sipms 
)

FUNCTIONS FOR SCIFICLASS ##################################################.

    It converts a particle position (an event) measured in cm to a position measured 
    in channels. The SiPM map is used. The position is in the scifi module frame.

Definition at line 276 of file ScifiCluster.py.

277 reverse=False, ch_max_num=ch_max_num, n_solid_sipms=n_solid_sipms):
278
279 """
280 It converts a particle position (an event) measured in cm to a position measured
281 in channels. The SiPM map is used. The position is in the scifi module frame.
282 """
283
284 if reverse is True:
285 for left, right, ch_range in sipm_map:
286 if left <= locpos <= right:
287 if not ch_range is False:
288 ch_start = ch_range[0]
289 return int(round(ch_max_num - ((locpos-left) / pitch + ch_start)))
290 else:
291 mapindex = int(round(n_solid_sipms - (locpos + hw) / charr))
292 return gaps_map.get(mapindex, False)
293
294 elif reverse is False:
295 for left, right, ch_range in sipm_map:
296 if left <= locpos <= right:
297 if not ch_range is False:
298 ch_start = ch_range[0]
299 return int(round(locpos-left) / pitch + ch_start)
300 else:
301 mapindex2 = int(round(locpos + hw) / charr)
302 return gaps_map.get( mapindex2, False)
303

◆ create_cluster()

ScifiCluster.create_cluster (   amplitude,
  width,
  wmp 
)
    The final function for creating a signal cluster

Definition at line 568 of file ScifiCluster.py.

568def create_cluster(amplitude, width, wmp):
569
570 """
571 The final function for creating a signal cluster
572
573 """
574 if wmp is False:
575 return False
576 if wmp is None:
577 return False
578
579 if not ((chpos_min < wmp) and (wmp< chpos_max)):
580 return False
581
582 if not ly_min <= amplitude >= width: return False
583 shifted_wmp = wmp + 0.5 # For right counting
584 cluster = cluster_generator(amplitude, width, shifted_wmp)
585
586 # Generate again if it doesn't look like real cluster
587 while is_realistic(cluster, width) is False:
588 cluster = cluster_generator(amplitude, width, shifted_wmp)
589
590 return cluster
591

◆ edep_to_ly()

ScifiCluster.edep_to_ly (   energy,
  CDF_integral = CDF_integral,
  energy_range = energy_range,
  ly_linear_params = ly_linear_params,
  k_cdfs_corr = k_cdfs_corr,
  sigma_in_percent = sigma_in_percent,
  ly_CDF_params = ly_CDF_params,
  CDF_ly_params = CDF_ly_params,
  ly_CDF_landau_params = ly_CDF_landau_params 
)
    It returns the light yield calculated from the energy deposit. The calculations are based 
    on the equality of the cumulative distribution functions (CDF) :
    energy => CDF(energy) => CDF(light yield) => ly

    The linear converting range 0.18 MeV < dE < 0.477 MeV corresponds 4.5 < LY < 104 ph.e.

    If energy more then 0.477 MeV the light yield calculated randomly (uniformly in the range) 
    according to the distribution

    Also a little randomness is added to the CDF value with a normal distribution and 
    standard deviation with 'sigma_in_percent' (in percent of the whole range 0 - max CDF)

Definition at line 478 of file ScifiCluster.py.

480 ly_CDF_params=ly_CDF_params, CDF_ly_params=CDF_ly_params, ly_CDF_landau_params=ly_CDF_landau_params):
481
482 """
483 It returns the light yield calculated from the energy deposit. The calculations are based
484 on the equality of the cumulative distribution functions (CDF) :
485 energy => CDF(energy) => CDF(light yield) => ly
486
487 The linear converting range 0.18 MeV < dE < 0.477 MeV corresponds 4.5 < LY < 104 ph.e.
488
489 If energy more then 0.477 MeV the light yield calculated randomly (uniformly in the range)
490 according to the distribution
491
492 Also a little randomness is added to the CDF value with a normal distribution and
493 standard deviation with 'sigma_in_percent' (in percent of the whole range 0 - max CDF)
494 """
495
496 e_min, e_max = energy_range
497 A, C = ly_linear_params
498 if e_min < energy < e_max:
499 ly_lin = A * energy + C
500 cdf_raw = approx_function(ly_lin, ly_CDF_landau_params) * k_cdfs_corr
501 elif e_max <= energy:
502 cdf_raw = CDF_integral * np.random.uniform(0., 1.0)
503 else:
504 return 0.
505 cdf_random = random.gauss(cdf_raw, sigma_in_percent * CDF_integral)
506
507 # Generate again while the CDF value is not in the range
508 while cdf_random < 0 or cdf_random > CDF_integral:
509 cdf_random = random.gauss(cdf_raw, sigma_in_percent * CDF_integral)
510
511 return approx_function(cdf_random, CDF_ly_params)
512

◆ fcheck_wall()

ScifiCluster.fcheck_wall (   cluster,
  wmp 
)
    A function that prevents clusters from propagating to the gaps between the mats

Definition at line 592 of file ScifiCluster.py.

592def fcheck_wall(cluster, wmp):
593
594 """
595 A function that prevents clusters from propagating to the gaps between the mats
596 """
597
598 if (wmp + 0.5)<= firstmatgap_startch: #1st mat
599 for chann in cluster:
600 if cluster[(chann + 0.5) > firstmatgap_startch] != 0.: cluster[chann] = 0.
601 return cluster
602 elif firstmatgap_endch <= (wmp + 0.5)<= secmatgap_startch: #2nd mat
603 for chann in cluster:
604 if (cluster[(chann + 0.5) < firstmatgap_endch] != 0.) or (cluster[(chann + 0.5) > secmatgap_startch] != 0.): cluster[chann] = 0.
605 return cluster
606 elif (wmp + 0.5) >= secmatgap_endch: #3rd mat
607 for chann in cluster:
608 if cluster[(chann + 0.5)< secmatgap_endch] != 0.: cluster[chann] = 0.
609 return cluster
610 else:
611 return False #"False" is returned if wmp lies in the gaps
612

◆ GetMatLength()

ScifiCluster.GetMatLength ( )
    It returns a length of a scifi mat. Identical for vertical and horizontal.

Definition at line 343 of file ScifiCluster.py.

343def GetMatLength():
344
345 """
346 It returns a length of a scifi mat. Identical for vertical and horizontal.
347 """
348
349 return scifimat_length
350

◆ GetMatNum()

ScifiCluster.GetMatNum (   DetID)
    It returns an id (number) of a scifi module. In current version one plane has 3 vertical
    and 3 horizontal scifi assemblies. 

Definition at line 334 of file ScifiCluster.py.

334def GetMatNum(DetID):
335
336 """
337 It returns an id (number) of a scifi module. In current version one plane has 3 vertical
338 and 3 horizontal scifi assemblies.
339 """
340
341 return int(DetID % 100 % 10)
342

◆ GetMatQty()

ScifiCluster.GetMatQty ( )
    It returns a number of scifi mats in a plane. Identical for vertical and horizontal fiber planes.

Definition at line 351 of file ScifiCluster.py.

351def GetMatQty():
352
353 """
354 It returns a number of scifi mats in a plane. Identical for vertical and horizontal fiber planes.
355 """
356
357 return n_planes
358

◆ GetMatType()

ScifiCluster.GetMatType (   DetID)
    It returns a type of a scifi module.
    1 - vertical scifi assembly
    0 - horizontal scifi assembly

Definition at line 324 of file ScifiCluster.py.

324def GetMatType(DetID):
325
326 """
327 It returns a type of a scifi module.
328 1 - vertical scifi assembly
329 0 - horizontal scifi assembly
330 """
331 return int((DetID % 100) // 10)
332
333

◆ GetSiPMs()

ScifiCluster.GetSiPMs ( )
    It returns a number of SiPMs per plane. Identical for vertical and horizontal fiber planes.

Definition at line 359 of file ScifiCluster.py.

359def GetSiPMs():
360
361 """
362 It returns a number of SiPMs per plane. Identical for vertical and horizontal fiber planes.
363 """
364
365 return n_sipms
366

◆ GetStationNum()

ScifiCluster.GetStationNum (   DetID)
    It returns an id of a plane. In current version the detector has 5 TT stations (numbering starts
    at 1.

Definition at line 367 of file ScifiCluster.py.

367def GetStationNum(DetID):
368
369 """
370 It returns an id of a plane. In current version the detector has 5 TT stations (numbering starts
371 at 1.
372 """
373 return int(DetID // 100)
374

◆ global_to_local()

ScifiCluster.global_to_local (   DetID,
  globalpos 
)

These functions assume that the detector is aligned with the beam ######.

    It returns the local coordinates in one scifi assembly frame from global coordinates.

Definition at line 379 of file ScifiCluster.py.

379def global_to_local(DetID, globalpos):
380
381 """
382 It returns the local coordinates in one scifi assembly frame from global coordinates.
383 """
384
385 global xshift, yshift
386
387 if GetMatType(DetID) == 0:
388 return globalpos - yshift
389
390 if GetMatType(DetID) == 1:
391 return globalpos - xshift
392
393

◆ is_realistic()

ScifiCluster.is_realistic (   cluster,
  width 
)
    It returns TRUE if cluster is realistic: it doesn't have a gap between numders, like
    [..., 0, 1, 2, 0, 0, 5, 6, ...], and it doens't have the light yield less then width.

Definition at line 553 of file ScifiCluster.py.

553def is_realistic(cluster, width):
554
555 """
556 It returns TRUE if cluster is realistic: it doesn't have a gap between numders, like
557 [..., 0, 1, 2, 0, 0, 5, 6, ...], and it doens't have the light yield less then width.
558 """
559
560 cluster_only_values = [(channel, value) for channel, value in enumerate(cluster) if value > 0]
561 first_channel, _ = cluster_only_values[0]
562 last_channel, _ = cluster_only_values[-1]
563 if len(cluster_only_values) != width or (last_channel - first_channel + 1) != width:
564 return False
565 else:
566 return True
567

◆ local_to_global()

ScifiCluster.local_to_global (   DetID,
  localpos 
)
    It returns the global coordinates from the local coordinates in one scifi assembly frame.

Definition at line 394 of file ScifiCluster.py.

394def local_to_global(DetID, localpos):
395
396 """
397 It returns the global coordinates from the local coordinates in one scifi assembly frame.
398 """
399
400 global xshift, yshift
401
402 if GetMatType(DetID) == 0:
403 return localpos + yshift
404
405 if GetMatType(DetID) == 1:
406 return localpos + xshift
407

◆ ly_attenuation()

ScifiCluster.ly_attenuation (   distance)
    It return the light yield losses in percent depending on the distance to SiPMs

Definition at line 419 of file ScifiCluster.py.

419def ly_attenuation(distance):
420
421 """
422 It return the light yield losses in percent depending on the distance to SiPMs
423 """
424
425 res_ly_loss = ly_loss_mean(distance) / ly_loss_mean(0.)
426 return res_ly_loss
427

◆ ly_loss_mean()

ScifiCluster.ly_loss_mean (   distance,
  params = ly_loss_params 
)
    It returns the light yield depending on the distance to SiPMs

Definition at line 410 of file ScifiCluster.py.

410def ly_loss_mean(distance, params=ly_loss_params):
411
412 """
413 It returns the light yield depending on the distance to SiPMs
414 """
415
416 A1, k1, A2, k2 = params
417 return A1 * exp(k1 * distance / 100.) + A2 * exp(k2 * distance / 100.)
418

◆ weigthed_mean_pos()

ScifiCluster.weigthed_mean_pos (   cluster)
    Calculate the weighted mean position of the cluster

Definition at line 613 of file ScifiCluster.py.

613def weigthed_mean_pos(cluster):
614
615 """
616 Calculate the weighted mean position of the cluster
617 """
618
619 if cluster is False:
620 return False
621 sumup = 0.
622 sumdown = 0.
623 wmp = 0.
624 for channel, value in enumerate(cluster):
625 if value > 0:
626 sumup += value * channel
627 sumdown += value
628 if sumdown != 0:
629 wmp = sumup / sumdown - 0.5
630 return wmp
631

Variable Documentation

◆ array

float ScifiCluster.array = gap + 2.*charr

Definition at line 43 of file ScifiCluster.py.

◆ biggap

float ScifiCluster.biggap = 0.038

Definition at line 44 of file ScifiCluster.py.

◆ CDF_integral

float ScifiCluster.CDF_integral = 0.0185640424

Definition at line 171 of file ScifiCluster.py.

◆ CDF_ly_params

dict ScifiCluster.CDF_ly_params
Initial value:
1= {
2 (0., 0.0006): (
3 (89., 4.152, 0.0001574, -1.326e+04, 4.3),
4 lambda cdf, *p: p[0] * sqrt(p[1]*(cdf+p[2])) * (1-exp(p[3]*cdf)) + p[4]
5 ),
6 (0.0006, 0.012): (
7 (158, 1.035, 0.24, 217),
8 lambda cdf, *p: p[0] * log(sqrt(cdf)+p[1]) + p[2]*exp(p[3]*cdf)
9 ),
10 (0.012, 0.018561405): (
11 (9.36, 335.984, -18100, -400, 15),
12 lambda cdf, *p: p[0] * log((p[1]-p[2]*cdf)/(p[1]+p[2]*cdf)) + p[3]*cdf + p[4]
13 ),
14 (0.018561405, 0.0185640424): (
15 (9.36, 335.984, -18100, -400, 15),
16 lambda cdf, *p: (p[0] * log((p[1]-p[2]*0.018561405)/(p[1]+p[2]*0.018561405))
17 + p[3]*0.018561405 + p[4])
18 )
19}

Definition at line 222 of file ScifiCluster.py.

◆ ch_max_num

int ScifiCluster.ch_max_num = 128*12 - 1

Definition at line 45 of file ScifiCluster.py.

◆ charr

float ScifiCluster.charr = 1.6

Definition at line 40 of file ScifiCluster.py.

◆ chpos_max

float ScifiCluster.chpos_max = 1535.5

Definition at line 145 of file ScifiCluster.py.

◆ chpos_min

float ScifiCluster.chpos_min = -0.5

Definition at line 144 of file ScifiCluster.py.

◆ cluster_width_max

int ScifiCluster.cluster_width_max = 10

Definition at line 147 of file ScifiCluster.py.

◆ cluster_width_mean_params

float ScifiCluster.cluster_width_mean_params = 0.6661, -0.006955, 2.163

Definition at line 155 of file ScifiCluster.py.

◆ cluster_width_min

int ScifiCluster.cluster_width_min = 1

Definition at line 146 of file ScifiCluster.py.

◆ cluster_width_sigma_params

float ScifiCluster.cluster_width_sigma_params = 1.103, -0.0005751

Definition at line 156 of file ScifiCluster.py.

◆ edge

float ScifiCluster.edge = 0.016

Definition at line 41 of file ScifiCluster.py.

◆ energy_range

float ScifiCluster.energy_range = 0.18, 0.477

Definition at line 168 of file ScifiCluster.py.

◆ firstmatgap_endch

float ScifiCluster.firstmatgap_endch = 523.5

Definition at line 139 of file ScifiCluster.py.

◆ firstmatgap_startch

float ScifiCluster.firstmatgap_startch = 521.5

Definition at line 138 of file ScifiCluster.py.

◆ gap

float ScifiCluster.gap = 0.022

Definition at line 42 of file ScifiCluster.py.

◆ gaps_map

dict ScifiCluster.gaps_map
Initial value:
1= {
2 0.: -0.5,
3 1.: 63.5,
4 2.: 127.5,
5 3.: 191.5,
6 4.: 255.5,
7 5.: 319.5,
8 6.: 383.5,
9 7.: 447.5,
10
11 8.: 511.5,
12 9.: 575.5,
13 10.: 639.5,
14 11.: 703.5,
15 12.: 767.5,
16 13.: 831.5,
17 14.: 895.5,
18 15.: 959.5,
19
20 16.: 1023.5,
21 17.: 1087.5,
22 18.: 1151.5,
23 19.: 1215.5,
24 20.: 1279.5,
25 21.: 1343.5,
26 22.: 1407.5,
27 23.: 1471.5,
28 24.: 1535.5
29}

Definition at line 105 of file ScifiCluster.py.

◆ hw

float ScifiCluster.hw = 19.528

Definition at line 38 of file ScifiCluster.py.

◆ is_sipm_hor_reversed

bool ScifiCluster.is_sipm_hor_reversed = True

Definition at line 23 of file ScifiCluster.py.

◆ is_sipm_vert_reversed

bool ScifiCluster.is_sipm_vert_reversed = True

Definition at line 24 of file ScifiCluster.py.

◆ k_cdfs_corr

float ScifiCluster.k_cdfs_corr = 0.993076766938

Definition at line 178 of file ScifiCluster.py.

◆ ly_CDF_landau_params

dict ScifiCluster.ly_CDF_landau_params
Initial value:
1= {
2 (4.5, 15): (
3 (0.001038, -0.000378, 3.53e-05),
4 lambda ly, *p: p[0] + p[1]*ly + p[2]*ly**2
5 ),
6 (15, 40): (
7 (-0.001986, -0.0003014, 7.031e-05, -2.067e-06, 1.892e-08),
8 lambda ly, *p: p[0] + p[1]*ly + p[2]*ly**2 + p[3]*ly**3 + p[4]*ly**4
9 ),
10 (40, 104): (
11 (-0.007149, 0.001056, -1.779e-05, 1.41e-07, -4.29e-10),
12 lambda ly, *p: p[0] + p[1]*ly + p[2]*ly**2 + p[3]*ly**3 + p[4]*ly**4
13 )
14}

Definition at line 206 of file ScifiCluster.py.

◆ ly_CDF_params

dict ScifiCluster.ly_CDF_params
Initial value:
1= {
2 (4.5, 13): (
3 (-13.2, 1.976),
4 lambda ly, *p: exp(p[0] + p[1]*sqrt(ly))
5 ),
6 (13, 104): (
7 (0.0108183, -0.179752, -19.2, 0.00772965),
8 lambda ly, *p: p[0]*(1 - exp(p[1]*(ly+p[2]))) / (1 + exp(p[1]*(ly+p[2]))) + p[3]
9 )
10}

Definition at line 192 of file ScifiCluster.py.

◆ ly_linear_params

float ScifiCluster.ly_linear_params = 332.882, -58.7085

Definition at line 175 of file ScifiCluster.py.

◆ ly_loss_params

float ScifiCluster.ly_loss_params = 20.78, -0.26, 7.89, -3.06

Definition at line 153 of file ScifiCluster.py.

◆ ly_loss_sigma_params

float ScifiCluster.ly_loss_sigma_params = 6.8482, -0.5757, 3.4458

Definition at line 154 of file ScifiCluster.py.

◆ ly_max

int ScifiCluster.ly_max = 104

Definition at line 164 of file ScifiCluster.py.

◆ ly_min

float ScifiCluster.ly_min = 4.5

Definition at line 163 of file ScifiCluster.py.

◆ n_planes

ScifiCluster.n_planes = snd_geo.Scifi.nscifi

Definition at line 259 of file ScifiCluster.py.

◆ n_sipms

ScifiCluster.n_sipms = snd_geo.Scifi.nsipms

Definition at line 260 of file ScifiCluster.py.

◆ n_solid_sipms

int ScifiCluster.n_solid_sipms = 12

Definition at line 46 of file ScifiCluster.py.

◆ nmats

ScifiCluster.nmats = snd_geo.Scifi.nmats

Definition at line 258 of file ScifiCluster.py.

◆ pitch

float ScifiCluster.pitch = 0.025

Definition at line 39 of file ScifiCluster.py.

◆ random_width_percent

float ScifiCluster.random_width_percent = 0.01

Definition at line 160 of file ScifiCluster.py.

◆ scifi_x

ScifiCluster.scifi_x = snd_geo.Scifi.xdim

Definition at line 254 of file ScifiCluster.py.

◆ scifi_y

ScifiCluster.scifi_y = snd_geo.Scifi.ydim

Definition at line 255 of file ScifiCluster.py.

◆ scifimat_gap

ScifiCluster.scifimat_gap = snd_geo.Scifi.scifimat_gap

Definition at line 252 of file ScifiCluster.py.

◆ scifimat_length

ScifiCluster.scifimat_length = snd_geo.Scifi.scifimat_length

Definition at line 251 of file ScifiCluster.py.

◆ scifimat_width

ScifiCluster.scifimat_width = snd_geo.Scifi.scifimat_width

Definition at line 250 of file ScifiCluster.py.

◆ scifimat_z

ScifiCluster.scifimat_z = snd_geo.Scifi.scifimat_z

Definition at line 256 of file ScifiCluster.py.

◆ secmatgap_endch

float ScifiCluster.secmatgap_endch = 1048.5

Definition at line 141 of file ScifiCluster.py.

◆ secmatgap_startch

float ScifiCluster.secmatgap_startch = 1046.5

Definition at line 140 of file ScifiCluster.py.

◆ sigma_from_width

int ScifiCluster.sigma_from_width = 1 / 4.

Definition at line 185 of file ScifiCluster.py.

◆ sigma_in_percent

float ScifiCluster.sigma_in_percent = 0.01

Definition at line 182 of file ScifiCluster.py.

◆ sipm_hor_pos

int ScifiCluster.sipm_hor_pos = +1

GLOBAL VALUES AND PARAMETERS ############.

Definition at line 18 of file ScifiCluster.py.

◆ sipm_map

tuple ScifiCluster.sipm_map

Definition at line 50 of file ScifiCluster.py.

◆ sipm_vert_pos

int ScifiCluster.sipm_vert_pos = -1

Definition at line 19 of file ScifiCluster.py.

◆ snd_geo

ScifiCluster.snd_geo = ConfigRegistry.loadpy("$SNDSW_ROOT/geometry/shipLHC_geom_config.py")

LOAD THE SCIFI GEOMETRY ############.

Definition at line 248 of file ScifiCluster.py.

◆ xshift

float ScifiCluster.xshift = -27.8

Need to automatize this part so it works for every geometry #########.

Definition at line 269 of file ScifiCluster.py.

◆ yshift

float ScifiCluster.yshift = +35.3

Definition at line 270 of file ScifiCluster.py.