Belle II Software development
TOPGeometryPar Class Reference

Singleton class for TOP Geometry Parameters. More...

#include <TOPGeometryPar.h>

Public Member Functions

virtual ~TOPGeometryPar ()
 Destructor.
 
void Initialize (const GearDir &content)
 Initialize from Gearbox (XML)
 
void Initialize ()
 Initialize from database.
 
bool isValid () const
 check if the geometry is available
 
const TOPGeometrygetGeometry () const
 Returns pointer to geometry object using basf2 units.
 
const FrontEndMappergetFrontEndMapper () const
 Returns front-end mapper (mapping of SCROD's to positions within TOP modules)
 
const ChannelMappergetChannelMapper () const
 Returns default channel mapper (mapping of channels to pixels)
 
const ChannelMappergetChannelMapper (ChannelMapper::EType type) const
 Returns channel mapper (mapping of channels to pixels) - Gearbox only.
 
double getPMTEfficiencyEnvelope (double energy) const
 Returns PMT efficiency envelope, e.g.
 
double getPMTEfficiency (double energy, int moduleID, int pmtID, double x, double y) const
 Returns PMT pixel efficiency, a product of quantum and collection efficiency.
 
double getRelativePixelEfficiency (int moduleID, int pixelID) const
 Returns relative pixel efficiency (including CE, RQE and threshold efficiency)
 
unsigned getPMTType (int moduleID, int pmtID) const
 Returns PMT type at a given position.
 
const TOPNominalTTSgetTTS (int moduleID, int pmtID) const
 Returns TTS of a PMT at given position.
 
double getPhaseIndex (double energy) const
 Returns phase refractive index of quartz at given photon energy.
 
double getGroupIndex (double energy) const
 Returns group refractive index of quartz at given photon energy.
 
double getPhaseIndexDerivative (double energy) const
 Returns the derivative (dn/dE) of phase refractive index of quartz at given photon energy.
 
double getGroupIndexDerivative (double energy) const
 Returns the derivative (dn_g/dE) of group refractive index of quartz at given photon energy.
 
double getAbsorptionLength (double energy) const
 Returns bulk absorption lenght of quartz at given photon energy.
 

Static Public Member Functions

static TOPGeometryParInstance ()
 Static method to obtain the pointer to its instance.
 

Static Public Attributes

static const double c_hc = 1239.84193
 Planck constant times speed of light in [eV*nm].
 

Private Member Functions

 TOPGeometryPar ()
 Hidden constructor since it is a singleton class.
 
void finalizeInitialization ()
 finalize initialization
 
TOPGeometrycreateConfiguration (const GearDir &content)
 Create a parameter object from gearbox.
 
TOPGeoBarSegment createBarSegment (const GearDir &content, const std::string &serialNumber)
 Create a parameter object from gearbox for bar segment.
 
TOPGeoMirrorSegment createMirrorSegment (const GearDir &content, const std::string &serialNumber)
 Create a parameter object from gearbox for mirror segment.
 
TOPGeoPrism createPrism (const GearDir &content, const std::string &serialNumber)
 Create a parameter object from gearbox for prism.
 
std::string addNumber (const std::string &str, unsigned number)
 Adds number to string.
 
void clearCache ()
 Clears cache for PMT dependent QE data - function is used in call backs.
 
void setEnvelopeQE () const
 Constructs envelope of quantum efficiency from PMT data.
 
void mapPmtQEToPositions () const
 Maps PMT QE data to positions within the detector.
 
void mapPmtTypeToPositions () const
 Maps PMT type to positions within the detector.
 
void prepareRelEfficiencies () const
 Prepares a map of relative pixel efficiencies.
 
int getUniquePmtID (int moduleID, int pmtID) const
 Returns unique PMT ID within the detector.
 
int getUniquePixelID (int moduleID, int pixelID) const
 Returns unique pixel ID within the detector.
 
double integralOfQE (const std::vector< float > &qe, double ce, double lambdaFirst, double lambdaStep) const
 Returns integral of quantum efficiency over photon energies.
 
double refractiveIndex (double lambda) const
 Quartz refractive index (SellMeier equation)
 

Private Attributes

TOPGeometrym_geo = 0
 geometry parameters from Gearbox
 
DBObjPtr< TOPGeometry > * m_geoDB = 0
 geometry parameters from database
 
bool m_fromDB = false
 parameters from database or Gearbox
 
bool m_valid = false
 true if geometry is available
 
bool m_oldPayload = false
 true if old payload found in DB
 
bool m_BfieldOn = true
 true if B field is on
 
FrontEndMapper m_frontEndMapper
 front end electronics mapper
 
ChannelMapper m_channelMapperIRS3B
 channel-pixel mapper
 
ChannelMapper m_channelMapperIRSX
 channel-pixel mapper
 
OptionalDBArray< TOPPmtInstallationm_pmtInstalled
 PMT installation data.
 
OptionalDBArray< TOPPmtQEm_pmtQEData
 quantum efficiencies
 
DBObjPtr< TOPCalChannelRQEm_channelRQE
 channel relative quantum effi.
 
DBObjPtr< TOPCalChannelThresholdEffm_thresholdEff
 channel threshold effi.
 
TOPNominalQE m_envelopeQE
 envelope quantum efficiency
 
std::map< int, const TOPPmtQE * > m_pmts
 QE data mapped to positions.
 
std::map< int, double > m_relEfficiencies
 pixel relative QE
 
std::map< int, unsigned > m_pmtTypes
 PMT types mapped to positions.
 

Static Private Attributes

static TOPGeometryPars_instance = 0
 Pointer to the class instance.
 

Detailed Description

Singleton class for TOP Geometry Parameters.

Definition at line 35 of file TOPGeometryPar.h.

Constructor & Destructor Documentation

◆ ~TOPGeometryPar()

~TOPGeometryPar ( )
virtual

Destructor.

Definition at line 35 of file TOPGeometryPar.cc.

36 {
37 if (m_geo) delete m_geo;
38 if (m_geoDB) delete m_geoDB;
39 s_instance = 0;
40 }
DBObjPtr< TOPGeometry > * m_geoDB
geometry parameters from database
TOPGeometry * m_geo
geometry parameters from Gearbox
static TOPGeometryPar * s_instance
Pointer to the class instance.

◆ TOPGeometryPar()

TOPGeometryPar ( )
inlineprivate

Hidden constructor since it is a singleton class.

Definition at line 182 of file TOPGeometryPar.h.

183 {}

Member Function Documentation

◆ addNumber()

std::string addNumber ( const std::string &  str,
unsigned  number 
)
private

Adds number to string.

Parameters
strstring
numbernumber to be added
Returns
string with a number

Definition at line 946 of file TOPGeometryPar.cc.

947 {
948 stringstream ss;
949 if (number < 10) {
950 ss << str << "0" << number;
951 } else {
952 ss << str << number;
953 }
954 string out;
955 ss >> out;
956 return out;
957 }

◆ clearCache()

void clearCache ( )
private

Clears cache for PMT dependent QE data - function is used in call backs.

Definition at line 157 of file TOPGeometryPar.cc.

158 {
160 m_pmts.clear();
161 m_relEfficiencies.clear();
162 m_pmtTypes.clear();
163 }
void clear()
Clears the object.
Definition: TOPNominalQE.h:75
TOPNominalQE m_envelopeQE
envelope quantum efficiency
std::map< int, unsigned > m_pmtTypes
PMT types mapped to positions.
std::map< int, const TOPPmtQE * > m_pmts
QE data mapped to positions.
std::map< int, double > m_relEfficiencies
pixel relative QE

◆ createBarSegment()

TOPGeoBarSegment createBarSegment ( const GearDir content,
const std::string &  serialNumber 
)
private

Create a parameter object from gearbox for bar segment.

Parameters
contentXML data directory
serialNumberbar segment serial number

Definition at line 868 of file TOPGeometryPar.cc.

870 {
871 // dimensions and material
872 GearDir params(content, "QuartzBars/QuartzBar[@SerialNumber='" + SN + "']");
873 TOPGeoBarSegment bar(params.getLength("Width"),
874 params.getLength("Thickness"),
875 params.getLength("Length"),
876 params.getString("Material"));
877 bar.setVendorData(params.getString("Vendor"), SN);
878
879 // optical surface
880 std::string surfaceName = params.getString("OpticalSurface");
881 double sigmaAlpha = params.getDouble("SigmaAlpha");
882 GearDir surfaceParams(content, "Modules/Surface[@name='" + surfaceName + "']");
883 auto& materials = geometry::Materials::getInstance();
884 auto quartzSurface = materials.createOpticalSurfaceConfig(surfaceParams);
885 bar.setSurface(quartzSurface, sigmaAlpha);
886
887 return bar;
888 }
static Materials & getInstance()
Get a reference to the singleton instance.
Definition: Materials.cc:85

◆ createConfiguration()

TOPGeometry * createConfiguration ( const GearDir content)
private

Create a parameter object from gearbox.

Parameters
contentXML data directory

Definition at line 437 of file TOPGeometryPar.cc.

438 {
439 TOPGeometry* geo = new TOPGeometry("TOPGeometry");
440
441 // PMT array
442
443 GearDir pmtParams(content, "PMTs/PMT");
444 TOPGeoPMT pmt(pmtParams.getLength("ModuleXSize"),
445 pmtParams.getLength("ModuleYSize"),
446 pmtParams.getLength("ModuleZSize") +
447 pmtParams.getLength("WindowThickness") +
448 pmtParams.getLength("BottomThickness"));
449 pmt.setWallThickness(pmtParams.getLength("ModuleWall"));
450 pmt.setWallMaterial(pmtParams.getString("wallMaterial"));
451 pmt.setFillMaterial(pmtParams.getString("fillMaterial"));
452 pmt.setSensVolume(pmtParams.getLength("SensXSize"),
453 pmtParams.getLength("SensYSize"),
454 pmtParams.getLength("SensThickness"),
455 pmtParams.getString("sensMaterial"));
456 pmt.setNumPixels(pmtParams.getInt("PadXNum"),
457 pmtParams.getInt("PadYNum"));
458 pmt.setWindow(pmtParams.getLength("WindowThickness"),
459 pmtParams.getString("winMaterial"));
460 pmt.setBottom(pmtParams.getLength("BottomThickness"),
461 pmtParams.getString("botMaterial"));
462
463 auto& materials = geometry::Materials::getInstance();
464 GearDir reflEdgeSurfParams(pmtParams, "reflectiveEdge/Surface");
465 pmt.setReflEdge(pmtParams.getLength("reflectiveEdge/width"),
466 pmtParams.getLength("reflectiveEdge/thickness"),
467 materials.createOpticalSurfaceConfig(reflEdgeSurfParams));
468
469 GearDir arrayParams(content, "PMTs");
470 TOPGeoPMTArray pmtArray(arrayParams.getInt("nPMTx"),
471 arrayParams.getInt("nPMTy"),
472 arrayParams.getLength("Xgap"),
473 arrayParams.getLength("Ygap"),
474 arrayParams.getString("stackMaterial"),
475 pmt);
476 pmtArray.setSiliconeCookie(arrayParams.getLength("siliconeCookie/thickness"),
477 arrayParams.getString("siliconeCookie/material"));
478 pmtArray.setWavelengthFilter(arrayParams.getLength("wavelengthFilter/thickness"),
479 arrayParams.getString("wavelengthFilter/material"));
480 pmtArray.setAirGap(arrayParams.getLength("airGap", 0));
481 double decoupledFraction = arrayParams.getDouble("decoupledFraction", 0);
482
483 // modules
484
485 GearDir moduleParams(content, "Modules");
486 GearDir glueParams(moduleParams, "Glue");
487 int numModules = moduleParams.getNumberNodes("Module");
488 for (int slotID = 1; slotID <= numModules; slotID++) {
489 std::string gearName = "Module[@slotID='" + std::to_string(slotID) + "']";
490 GearDir slotParams(moduleParams, gearName);
491 TOPGeoModule module(slotID,
492 slotParams.getLength("Radius"),
493 slotParams.getAngle("Phi"),
494 slotParams.getLength("BackwardZ"));
495 int cNumber = slotParams.getInt("ConstructionNumber");
496 module.setModuleCNumber(cNumber);
497 module.setName(addNumber(module.getName(), cNumber));
498
499 auto prism = createPrism(content, slotParams.getString("Prism"));
500 prism.setName(addNumber(prism.getName(), cNumber));
501 module.setPrism(prism);
502
503 auto barSegment2 = createBarSegment(content, slotParams.getString("BarSegment2"));
504 barSegment2.setName(addNumber(barSegment2.getName() + "2-", cNumber));
505 barSegment2.setGlue(glueParams.getLength("Thicknes1"),
506 glueParams.getString("GlueMaterial"));
507 module.setBarSegment2(barSegment2);
508
509 auto barSegment1 = createBarSegment(content, slotParams.getString("BarSegment1"));
510 barSegment1.setName(addNumber(barSegment1.getName() + "1-", cNumber));
511 barSegment1.setGlue(glueParams.getLength("Thicknes2"),
512 glueParams.getString("GlueMaterial"));
513 module.setBarSegment1(barSegment1);
514
515 auto mirror = createMirrorSegment(content, slotParams.getString("Mirror"));
516 mirror.setName(addNumber(mirror.getName(), cNumber));
517 mirror.setGlue(glueParams.getLength("Thicknes3"),
518 glueParams.getString("GlueMaterial"));
519 module.setMirrorSegment(mirror);
520
521 module.setPMTArray(pmtArray);
522 if (decoupledFraction > 0) module.generateDecoupledPMTs(decoupledFraction);
523
524 geo->appendModule(module);
525 }
526
527 // displaced geometry (if defined)
528
529 GearDir displacedGeometry(content, "DisplacedGeometry");
530 if (displacedGeometry) {
531 if (displacedGeometry.getInt("SwitchON") != 0) {
532 B2WARNING("TOP: displaced geometry is activated");
533 for (const GearDir& slot : displacedGeometry.getNodes("Slot")) {
534 int moduleID = slot.getInt("@ID");
535 if (!geo->isModuleIDValid(moduleID)) {
536 B2WARNING("TOPGeometryPar: DisplacedGeometry.xml: invalid moduleID."
537 << LogVar("moduleID", moduleID));
538 continue;
539 }
540 TOPGeoModuleDisplacement moduleDispl(slot.getLength("x"),
541 slot.getLength("y"),
542 slot.getLength("z"),
543 slot.getAngle("alpha"),
544 slot.getAngle("beta"),
545 slot.getAngle("gamma"));
546 auto& module = const_cast<TOPGeoModule&>(geo->getModule(moduleID));
547 module.setModuleDisplacement(moduleDispl);
548 }
549 }
550 }
551
552 // displaced PMT arrays (if defined)
553
554 GearDir displacedPMTArrays(content, "DisplacedPMTArrays");
555 if (displacedPMTArrays) {
556 if (displacedPMTArrays.getInt("SwitchON") != 0) {
557 B2WARNING("TOP: displaced PMT arrays are activated");
558 for (const GearDir& slot : displacedPMTArrays.getNodes("Slot")) {
559 int moduleID = slot.getInt("@ID");
560 if (!geo->isModuleIDValid(moduleID)) {
561 B2WARNING("TOPGeometryPar: DisplacedPMTArrays.xml: invalid moduleID."
562 << LogVar("moduleID", moduleID));
563 continue;
564 }
565 TOPGeoPMTArrayDisplacement arrayDispl(slot.getLength("x"),
566 slot.getLength("y"),
567 slot.getAngle("alpha"));
568 auto& module = const_cast<TOPGeoModule&>(geo->getModule(moduleID));
569 module.setPMTArrayDisplacement(arrayDispl);
570 }
571 }
572 }
573
574 // broken glues (if any)
575
576 GearDir brokenGlues(content, "BrokenGlues");
577 if (brokenGlues) {
578 if (brokenGlues.getInt("SwitchON") != 0) {
579 auto material = brokenGlues.getString("material");
580 for (const GearDir& slot : brokenGlues.getNodes("Slot")) {
581 int moduleID = slot.getInt("@ID");
582 if (!geo->isModuleIDValid(moduleID)) {
583 B2WARNING("TOPGeometryPar: BrokenGlues.xml: invalid moduleID."
584 << LogVar("moduleID", moduleID));
585 continue;
586 }
587 auto& module = const_cast<TOPGeoModule&>(geo->getModule(moduleID));
588 for (const GearDir& glue : slot.getNodes("Glue")) {
589 int glueID = glue.getInt("@ID");
590 double fraction = glue.getDouble("fraction");
591 if (fraction <= 0) continue;
592 double angle = glue.getAngle("angle");
593 module.setBrokenGlue(glueID, fraction, angle, material);
594 }
595 }
596 }
597 }
598
599 // peel-off cookies (if any)
600
601 GearDir peelOff(content, "PeelOffCookies");
602 if (peelOff) {
603 if (peelOff.getInt("SwitchON") != 0) {
604 auto material = peelOff.getString("material");
605 double thickness = peelOff.getLength("thickness");
606 for (const GearDir& slot : peelOff.getNodes("Slot")) {
607 int moduleID = slot.getInt("@ID");
608 if (!geo->isModuleIDValid(moduleID)) {
609 B2WARNING("TOPGeometryPar: PeelOffCookiess.xml: invalid moduleID."
610 << LogVar("moduleID", moduleID));
611 continue;
612 }
613 auto& module = const_cast<TOPGeoModule&>(geo->getModule(moduleID));
614 module.setPeelOffRegions(thickness, material);
615 for (const GearDir& region : slot.getNodes("Region")) {
616 int regionID = region.getInt("@ID");
617 double fraction = region.getDouble("fraction");
618 if (fraction <= 0) continue;
619 double angle = region.getAngle("angle");
620 module.appendPeelOffRegion(regionID, fraction, angle);
621 }
622 }
623 }
624 }
625
626 // front-end electronics geometry
627
628 GearDir feParams(content, "FrontEndGeo");
629 GearDir fbParams(feParams, "FrontBoard");
630 TOPGeoFrontEnd frontEnd;
631 frontEnd.setFrontBoard(fbParams.getLength("width"),
632 fbParams.getLength("height"),
633 fbParams.getLength("thickness"),
634 fbParams.getLength("gap"),
635 fbParams.getLength("y"),
636 fbParams.getString("material"));
637 GearDir hvParams(feParams, "HVBoard");
638 frontEnd.setHVBoard(hvParams.getLength("width"),
639 hvParams.getLength("length"),
640 hvParams.getLength("thickness"),
641 hvParams.getLength("gap"),
642 hvParams.getLength("y"),
643 hvParams.getString("material"));
644 GearDir bsParams(feParams, "BoardStack");
645 frontEnd.setBoardStack(bsParams.getLength("width"),
646 bsParams.getLength("height"),
647 bsParams.getLength("length"),
648 bsParams.getLength("gap"),
649 bsParams.getLength("y"),
650 bsParams.getString("material"),
651 bsParams.getLength("spacerWidth"),
652 bsParams.getString("spacerMaterial"));
653 geo->setFrontEnd(frontEnd, feParams.getInt("numBoardStacks"));
654
655 // QBB
656
657 GearDir qbbParams(content, "QBB");
658 TOPGeoQBB qbb(qbbParams.getLength("width"),
659 qbbParams.getLength("length"),
660 qbbParams.getLength("prismPosition"),
661 qbbParams.getString("material"));
662
663 GearDir outerPanelParams(qbbParams, "outerPanel");
664 TOPGeoHoneycombPanel outerPanel(outerPanelParams.getLength("width"),
665 outerPanelParams.getLength("length"),
666 outerPanelParams.getLength("minThickness"),
667 outerPanelParams.getLength("maxThickness"),
668 outerPanelParams.getLength("radius"),
669 outerPanelParams.getLength("edgeWidth"),
670 outerPanelParams.getLength("y"),
671 outerPanelParams.getInt("N"),
672 outerPanelParams.getString("material"),
673 outerPanelParams.getString("edgeMaterial"),
674 "TOPOuterHoneycombPanel");
675 qbb.setOuterPanel(outerPanel);
676
677 GearDir innerPanelParams(qbbParams, "innerPanel");
678 TOPGeoHoneycombPanel innerPanel(innerPanelParams.getLength("width"),
679 innerPanelParams.getLength("length"),
680 innerPanelParams.getLength("minThickness"),
681 innerPanelParams.getLength("maxThickness"),
682 innerPanelParams.getLength("radius"),
683 innerPanelParams.getLength("edgeWidth"),
684 innerPanelParams.getLength("y"),
685 innerPanelParams.getInt("N"),
686 innerPanelParams.getString("material"),
687 innerPanelParams.getString("edgeMaterial"),
688 "TOPInnerHoneycombPanel");
689 qbb.setInnerPanel(innerPanel);
690
691 GearDir sideRailsParams(qbbParams, "sideRails");
692 TOPGeoSideRails sideRails(sideRailsParams.getLength("thickness"),
693 sideRailsParams.getLength("reducedThickness"),
694 sideRailsParams.getLength("height"),
695 sideRailsParams.getString("material"));
696 qbb.setSideRails(sideRails);
697
698 GearDir prismEnclParams(qbbParams, "prismEnclosure");
699 TOPGeoPrismEnclosure prismEncl(prismEnclParams.getLength("length"),
700 prismEnclParams.getLength("height"),
701 prismEnclParams.getAngle("angle"),
702 prismEnclParams.getLength("bottomThickness"),
703 prismEnclParams.getLength("sideThickness"),
704 prismEnclParams.getLength("backThickness"),
705 prismEnclParams.getLength("frontThickness"),
706 prismEnclParams.getLength("extensionThickness"),
707 prismEnclParams.getString("material"));
708 qbb.setPrismEnclosure(prismEncl);
709
710 GearDir endPlateParams(qbbParams, "forwardEndPlate");
711 TOPGeoEndPlate endPlate(endPlateParams.getLength("thickness"),
712 endPlateParams.getLength("height"),
713 endPlateParams.getString("material"),
714 "TOPForwardEndPlate");
715 qbb.setEndPlate(endPlate);
716
717 GearDir coldPlateParams(qbbParams, "coldPlate");
718 TOPGeoColdPlate coldPlate(coldPlateParams.getLength("baseThickness"),
719 coldPlateParams.getString("baseMaterial"),
720 coldPlateParams.getLength("coolThickness"),
721 coldPlateParams.getLength("coolWidth"),
722 coldPlateParams.getString("coolMaterial"));
723 qbb.setColdPlate(coldPlate);
724
725 geo->setQBB(qbb);
726
727 // nominal QE
728
729 GearDir qeParams(content, "QE");
730 std::vector<float> qeData;
731 for (const GearDir& Qeffi : qeParams.getNodes("Qeffi")) {
732 qeData.push_back(Qeffi.getDouble(""));
733 }
734 TOPNominalQE nominalQE(qeParams.getLength("LambdaFirst") / Unit::nm,
735 qeParams.getLength("LambdaStep") / Unit::nm,
736 qeParams.getDouble("ColEffi"),
737 qeData);
738 geo->setNominalQE(nominalQE);
739
740 // nominal TTS
741
742 GearDir ttsParams(content, "PMTs/TTS");
743 TOPNominalTTS nominalTTS("TOPNominalTTS");
744 for (const GearDir& Gauss : ttsParams.getNodes("Gauss")) {
745 nominalTTS.appendGaussian(Gauss.getDouble("fraction"),
746 Gauss.getTime("mean"),
747 Gauss.getTime("sigma"));
748 }
749 nominalTTS.normalize();
750 geo->setNominalTTS(nominalTTS);
751
752 // PMT type dependent TTS
753
754 GearDir pmtTTSParams(content, "TTSofPMTs");
755 for (const GearDir& ttsPar : pmtTTSParams.getNodes("TTSpar")) {
756 int type = ttsPar.getInt("type");
757 double tuneFactor = ttsPar.getDouble("PDEtuneFactor");
758 TOPNominalTTS tts("TTS of " + ttsPar.getString("@name") + " PMT");
759 tts.setPMTType(type);
760 for (const GearDir& Gauss : ttsPar.getNodes("Gauss")) {
761 tts.appendGaussian(Gauss.getDouble("fraction"),
762 Gauss.getTime("mean"),
763 Gauss.getTime("sigma"));
764 }
765 tts.normalize();
766 geo->appendTTS(tts);
767 geo->appendPDETuningFactor(type, tuneFactor);
768 }
769
770 // nominal TDC
771
772 GearDir tdcParams(content, "TDC");
773 if (tdcParams) {
774 TOPNominalTDC nominalTDC(tdcParams.getInt("numWindows"),
775 tdcParams.getInt("subBits"),
776 tdcParams.getTime("syncTimeBase"),
777 tdcParams.getInt("numofBunches"),
778 tdcParams.getTime("offset"),
779 tdcParams.getTime("pileupTime"),
780 tdcParams.getTime("doubleHitResolution"),
781 tdcParams.getTime("timeJitter"),
782 tdcParams.getDouble("efficiency"));
783 nominalTDC.setADCBits(tdcParams.getInt("adcBits"));
784 nominalTDC.setAveragePedestal(tdcParams.getInt("averagePedestal"));
785 geo->setNominalTDC(nominalTDC);
786 } else {
787 TOPNominalTDC nominalTDC(pmtParams.getInt("TDCbits"),
788 pmtParams.getTime("TDCbitwidth"),
789 pmtParams.getTime("TDCoffset", 0),
790 pmtParams.getTime("TDCpileupTime", 0),
791 pmtParams.getTime("TDCdoubleHitResolution", 0),
792 pmtParams.getTime("TDCtimeJitter", 50e-3),
793 pmtParams.getDouble("TDCefficiency", 1));
794 geo->setNominalTDC(nominalTDC);
795 }
796
797 // single photon signal shape
798
799 GearDir shapeParams(content, "SignalShape");
800 if (shapeParams) {
801 GearDir noiseBandwidth(shapeParams, "noiseBandwidth");
802 TOPSignalShape signalShape(shapeParams.getArray("sampleValues"),
803 geo->getNominalTDC().getSampleWidth(),
804 shapeParams.getTime("tailTimeConstant"),
805 noiseBandwidth.getDouble("pole1") / 1000,
806 noiseBandwidth.getDouble("pole2") / 1000);
807 geo->setSignalShape(signalShape);
808 }
809
810 // calibration pulse shape
811
812 GearDir calpulseParams(content, "CalPulseShape");
813 if (calpulseParams) {
814 GearDir noiseBandwidth(calpulseParams, "noiseBandwidth");
815 TOPSignalShape shape(calpulseParams.getArray("sampleValues"),
816 geo->getNominalTDC().getSampleWidth(),
817 calpulseParams.getTime("tailTimeConstant"),
818 noiseBandwidth.getDouble("pole1") / 1000,
819 noiseBandwidth.getDouble("pole2") / 1000);
820 geo->setCalPulseShape(shape);
821 }
822
823 // wavelength filter bulk transmittance
824
825 std::string materialNode = "Materials/Material[@name='TOPWavelengthFilterIHU340']";
826 GearDir filterMaterial(content, materialNode);
827 if (!filterMaterial) {
828 B2FATAL("TOPGeometry: " << materialNode << " not found");
829 }
830 GearDir property(filterMaterial, "Property[@name='ABSLENGTH']");
831 if (!property) {
832 B2FATAL("TOPGeometry: " << materialNode << ", Property ABSLENGTH not found");
833 }
834 int numNodes = property.getNumberNodes("value");
835 if (numNodes > 1) {
836 double conversion = Unit::convertValue(1, property.getString("@unit", "GeV"));
837 std::vector<double> energies;
838 std::vector<double> absLengths;
839 for (int i = 0; i < numNodes; i++) {
840 GearDir value(property, "value", i + 1);
841 energies.push_back(value.getDouble("@energy") * conversion / Unit::eV);// [eV]
842 absLengths.push_back(value.getDouble() * Unit::mm); // [cm]
843 }
844 TSpline3 spline("absLen", energies.data(), absLengths.data(), energies.size());
845 double lambdaFirst = c_hc / energies.back();
846 double lambdaLast = c_hc / energies[0];
847 double lambdaStep = 5; // [nm]
848 int numSteps = (lambdaLast - lambdaFirst) / lambdaStep + 1;
849 const double filterThickness = arrayParams.getLength("wavelengthFilter/thickness");
850 std::vector<float> bulkTransmittances;
851 for (int i = 0; i < numSteps; i++) {
852 double wavelength = lambdaFirst + lambdaStep * i;
853 double energy = c_hc / wavelength;
854 double absLen = spline.Eval(energy);
855 bulkTransmittances.push_back(exp(-filterThickness / absLen));
856 }
857 TOPWavelengthFilter filter(lambdaFirst, lambdaStep, bulkTransmittances);
858 geo->setWavelengthFilter(filter);
859 } else {
860 B2FATAL("TOPGeometry: " << materialNode
861 << ", Property ABSLENGTH has less than 2 nodes");
862 }
863
864 return geo;
865 }
TOPGeoPrism createPrism(const GearDir &content, const std::string &serialNumber)
Create a parameter object from gearbox for prism.
TOPGeoMirrorSegment createMirrorSegment(const GearDir &content, const std::string &serialNumber)
Create a parameter object from gearbox for mirror segment.
std::string addNumber(const std::string &str, unsigned number)
Adds number to string.
static const double c_hc
Planck constant times speed of light in [eV*nm].
TOPGeoBarSegment createBarSegment(const GearDir &content, const std::string &serialNumber)
Create a parameter object from gearbox for bar segment.
static const double mm
[millimeters]
Definition: Unit.h:70
static const double nm
[nanometers]
Definition: Unit.h:72
static const double eV
[electronvolt]
Definition: Unit.h:112
Class to store variables with their name which were sent to the logging service.
static double convertValue(double value, const std::string &unitString)
Converts a floating point value to the standard framework unit.
Definition: UnitConst.cc:129
std::map< ExpRun, std::pair< double, double > > filter(const std::map< ExpRun, std::pair< double, double > > &runs, double cut, std::map< ExpRun, std::pair< double, double > > &runsRemoved)
filter events to remove runs shorter than cut, it stores removed runs in runsRemoved
Definition: Splitter.cc:38

◆ createMirrorSegment()

TOPGeoMirrorSegment createMirrorSegment ( const GearDir content,
const std::string &  serialNumber 
)
private

Create a parameter object from gearbox for mirror segment.

Parameters
contentXML data directory
serialNumbermirror segment serial number

Definition at line 891 of file TOPGeometryPar.cc.

893 {
894 // dimensions and material
895 GearDir params(content, "Mirrors/Mirror[@SerialNumber='" + SN + "']");
896 TOPGeoMirrorSegment mirror(params.getLength("Width"),
897 params.getLength("Thickness"),
898 params.getLength("Length"),
899 params.getString("Material"));
900 mirror.setVendorData(params.getString("Vendor"), SN);
901 mirror.setRadius(params.getLength("Radius"));
902 mirror.setCenterOfCurvature(params.getLength("Xpos"), params.getLength("Ypos"));
903
904 // mirror reflective coating
905 auto& materials = geometry::Materials::getInstance();
906 GearDir coatingParams(params, "Surface");
907 mirror.setCoating(params.getLength("mirrorThickness"), "Al",
908 materials.createOpticalSurfaceConfig(coatingParams));
909
910 // optical surface
911 std::string surfaceName = params.getString("OpticalSurface");
912 double sigmaAlpha = params.getDouble("SigmaAlpha");
913 GearDir surfaceParams(content, "Modules/Surface[@name='" + surfaceName + "']");
914 auto quartzSurface = materials.createOpticalSurfaceConfig(surfaceParams);
915 mirror.setSurface(quartzSurface, sigmaAlpha);
916
917 return mirror;
918 }

◆ createPrism()

TOPGeoPrism createPrism ( const GearDir content,
const std::string &  serialNumber 
)
private

Create a parameter object from gearbox for prism.

Parameters
contentXML data directory
serialNumberprism serial number

Definition at line 921 of file TOPGeometryPar.cc.

923 {
924 // dimensions and material
925 GearDir params(content, "Prisms/Prism[@SerialNumber='" + SN + "']");
926 TOPGeoPrism prism(params.getLength("Width"),
927 params.getLength("Thickness"),
928 params.getLength("Length"),
929 params.getLength("ExitThickness"),
930 0.,
931 params.getString("Material"));
932 prism.setAngle(params.getAngle("Angle"));
933 prism.setVendorData(params.getString("Vendor"), SN);
934
935 // optical surface
936 std::string surfaceName = params.getString("OpticalSurface");
937 double sigmaAlpha = params.getDouble("SigmaAlpha");
938 GearDir surfaceParams(content, "Modules/Surface[@name='" + surfaceName + "']");
939 auto& materials = geometry::Materials::getInstance();
940 auto quartzSurface = materials.createOpticalSurfaceConfig(surfaceParams);
941 prism.setSurface(quartzSurface, sigmaAlpha);
942
943 return prism;
944 }

◆ finalizeInitialization()

void finalizeInitialization ( )
private

finalize initialization

Definition at line 135 of file TOPGeometryPar.cc.

136 {
137 // set B field flag
138 m_BfieldOn = (BFieldManager::getField(0, 0, 0).R() / Unit::T) > 0.1;
139
140 // add call backs for PMT data
141 m_pmtInstalled.addCallback(this, &TOPGeometryPar::clearCache);
142 m_pmtQEData.addCallback(this, &TOPGeometryPar::clearCache);
143
144 // print geometry if the debug level for 'top' is set 10000
145 const auto& logSystem = LogSystem::Instance();
146 if (logSystem.isLevelEnabled(LogConfig::c_Debug, 10000, "top")) {
147 getGeometry()->print();
148 if (m_oldPayload) {
149 cout << "Envelope QE same as nominal quantum efficiency" << endl << endl;
150 return;
151 }
153 m_envelopeQE.print("Envelope QE");
154 }
155 }
@ c_Debug
Debug: for code development.
Definition: LogConfig.h:26
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
Definition: LogSystem.cc:31
bool isEmpty() const
Checks the status.
Definition: TOPNominalQE.h:88
OptionalDBArray< TOPPmtQE > m_pmtQEData
quantum efficiencies
bool m_BfieldOn
true if B field is on
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
void clearCache()
Clears cache for PMT dependent QE data - function is used in call backs.
bool m_oldPayload
true if old payload found in DB
void setEnvelopeQE() const
Constructs envelope of quantum efficiency from PMT data.
OptionalDBArray< TOPPmtInstallation > m_pmtInstalled
PMT installation data.
static const double T
[tesla]
Definition: Unit.h:120
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91
void print(const std::string &title="Nominal quantum efficiency") const override
Print the content of the class.
Definition: TOPNominalQE.cc:68
void print(const std::string &title="TOP geometry parameters") const override
Print the content of the class.
Definition: TOPGeometry.cc:142

◆ getAbsorptionLength()

double getAbsorptionLength ( double  energy) const
inline

Returns bulk absorption lenght of quartz at given photon energy.

Parameters
energyphoton energy [eV]
Returns
bulk absorption lenght

Definition at line 340 of file TOPGeometryPar.h.

341 {
342 double lambda = c_hc / energy;
343 return 15100 * pow(lambda / 405, 4); // Alan Schwartz, 2013 (private comunication)
344 }

◆ getChannelMapper() [1/2]

const ChannelMapper & getChannelMapper ( ) const
inline

Returns default channel mapper (mapping of channels to pixels)

Returns
channel mapper object

Definition at line 83 of file TOPGeometryPar.h.

83{return m_channelMapperIRSX;}
ChannelMapper m_channelMapperIRSX
channel-pixel mapper

◆ getChannelMapper() [2/2]

const ChannelMapper & getChannelMapper ( ChannelMapper::EType  type) const
inline

Returns channel mapper (mapping of channels to pixels) - Gearbox only.

Returns
channel mapper object

Definition at line 89 of file TOPGeometryPar.h.

90 {
91 switch (type) {
92 case ChannelMapper::c_IRS3B: return m_channelMapperIRS3B;
93 case ChannelMapper::c_IRSX: return m_channelMapperIRSX;
94 default: return m_channelMapperIRSX;
95 }
96 }
ChannelMapper m_channelMapperIRS3B
channel-pixel mapper

◆ getFrontEndMapper()

const FrontEndMapper & getFrontEndMapper ( ) const
inline

Returns front-end mapper (mapping of SCROD's to positions within TOP modules)

Returns
front-end mapper object

Definition at line 77 of file TOPGeometryPar.h.

77{return m_frontEndMapper;}
FrontEndMapper m_frontEndMapper
front end electronics mapper

◆ getGeometry()

const TOPGeometry * getGeometry ( ) const

Returns pointer to geometry object using basf2 units.

Returns
pointer to geometry object

Definition at line 165 of file TOPGeometryPar.cc.

166 {
167 if (!m_valid) B2FATAL("No geometry available for TOP");
168
170 if (m_fromDB) {
171 return &(**m_geoDB);
172 } else {
173 return m_geo;
174 }
175 }
static void useBasf2Units()
Use basf2 units when returning geometry parameters.
Definition: TOPGeometry.h:53
bool m_fromDB
parameters from database or Gearbox
bool m_valid
true if geometry is available

◆ getGroupIndex()

double getGroupIndex ( double  energy) const

Returns group refractive index of quartz at given photon energy.

Parameters
energyphoton energy [eV]
Returns
group refractive index

Definition at line 981 of file TOPGeometryPar.cc.

982 {
983 double lambda = c_hc / energy;
984 double dl = 1.0; // [nm]
985 double n = refractiveIndex(lambda);
986 double dndl = (refractiveIndex(lambda + dl / 2) - refractiveIndex(lambda - dl / 2)) / dl;
987 return n / (1 + lambda / n * dndl);
988 }
double refractiveIndex(double lambda) const
Quartz refractive index (SellMeier equation)

◆ getGroupIndexDerivative()

double getGroupIndexDerivative ( double  energy) const
inline

Returns the derivative (dn_g/dE) of group refractive index of quartz at given photon energy.

Parameters
energyphoton energy [eV]
Returns
group refractive index

Definition at line 334 of file TOPGeometryPar.h.

335 {
336 double dE = 0.01; // [eV]
337 return (getGroupIndex(energy + dE / 2) - getGroupIndex(energy - dE / 2)) / dE;
338 }
double getGroupIndex(double energy) const
Returns group refractive index of quartz at given photon energy.

◆ getPhaseIndex()

double getPhaseIndex ( double  energy) const

Returns phase refractive index of quartz at given photon energy.

Parameters
energyphoton energy [eV]
Returns
phase refractive index

Definition at line 975 of file TOPGeometryPar.cc.

976 {
977 double lambda = c_hc / energy;
978 return refractiveIndex(lambda);
979 }

◆ getPhaseIndexDerivative()

double getPhaseIndexDerivative ( double  energy) const
inline

Returns the derivative (dn/dE) of phase refractive index of quartz at given photon energy.

Parameters
energyphoton energy [eV]
Returns
derivative of phase refractive index

Definition at line 328 of file TOPGeometryPar.h.

329 {
330 double dE = 0.01; // [eV]
331 return (getPhaseIndex(energy + dE / 2) - getPhaseIndex(energy - dE / 2)) / dE;
332 }
double getPhaseIndex(double energy) const
Returns phase refractive index of quartz at given photon energy.

◆ getPMTEfficiency()

double getPMTEfficiency ( double  energy,
int  moduleID,
int  pmtID,
double  x,
double  y 
) const

Returns PMT pixel efficiency, a product of quantum and collection efficiency.

Parameters
energyphoton energy in [eV]
moduleIDslot ID
pmtIDPMT ID
xphoton detection position x in local PMT frame
yphoton detection position y in local PMT frame
Returns
the efficiency

Definition at line 189 of file TOPGeometryPar.cc.

192 {
193 const auto* geo = getGeometry();
194 if (!geo->isModuleIDValid(moduleID)) return 0;
195
196 double lambda = c_hc / energy;
197
198 if (m_oldPayload) { // filter transmittance is included in nominal QE, return it!
199 return geo->getNominalQE().getEfficiency(lambda);
200 }
201
202 if (m_pmts.empty()) mapPmtQEToPositions();
203
204 int id = getUniquePmtID(moduleID, pmtID);
205 const auto* pmtQE = m_pmts[id];
206 if (!pmtQE) return 0;
207
208 const auto& pmtArray = geo->getModule(moduleID).getPMTArray();
209 auto pmtPixel = pmtArray.getPMT().getPixelID(x, y);
210 if (pmtPixel == 0) return 0;
211
212 auto pixelID = pmtArray.getPixelID(pmtID, pmtPixel);
213 auto channel = getChannelMapper().getChannel(pixelID);
214
215 double RQE = geo->getPDETuningFactor(getPMTType(moduleID, pmtID));
216 if (m_channelRQE.isValid()) RQE *= m_channelRQE->getRQE(moduleID, channel);
217
218 return pmtQE->getEfficiency(pmtPixel, lambda, m_BfieldOn) * RQE;
219
220 }
unsigned getChannel(int pixel) const
Converts pixel to hardware channel number (0-based)
unsigned getPMTType(int moduleID, int pmtID) const
Returns PMT type at a given position.
const ChannelMapper & getChannelMapper() const
Returns default channel mapper (mapping of channels to pixels)
DBObjPtr< TOPCalChannelRQE > m_channelRQE
channel relative quantum effi.
void mapPmtQEToPositions() const
Maps PMT QE data to positions within the detector.
int getUniquePmtID(int moduleID, int pmtID) const
Returns unique PMT ID within the detector.

◆ getPMTEfficiencyEnvelope()

double getPMTEfficiencyEnvelope ( double  energy) const

Returns PMT efficiency envelope, e.g.

at given photon energy the maximum over all PMT's of a product of quantum and collection efficiency.

Parameters
energyphoton energy in [eV]
Returns
the maximal efficiency

Definition at line 177 of file TOPGeometryPar.cc.

178 {
179 double lambda = c_hc / energy;
180
181 if (m_oldPayload) { // filter transmittance is included in nominal QE, return it!
182 return getGeometry()->getNominalQE().getEfficiency(lambda);
183 }
184
186 return m_envelopeQE.getEfficiency(lambda);
187 }
const TOPNominalQE & getNominalQE() const
Returns nominal quantum efficiency of PMT.
Definition: TOPGeometry.h:180
double getEfficiency(double lambda) const
Returns quantum times collection efficiency at given photon wavelength using linear interpolation.
Definition: TOPNominalQE.h:109

◆ getPMTType()

unsigned getPMTType ( int  moduleID,
int  pmtID 
) const

Returns PMT type at a given position.

Parameters
moduleIDslot ID
pmtIDPMT ID
Returns
PMT type

Definition at line 246 of file TOPGeometryPar.cc.

247 {
248 if (m_pmtTypes.empty()) mapPmtTypeToPositions();
249
250 int id = getUniquePmtID(moduleID, pmtID);
251 return m_pmtTypes[id];
252 }
void mapPmtTypeToPositions() const
Maps PMT type to positions within the detector.

◆ getRelativePixelEfficiency()

double getRelativePixelEfficiency ( int  moduleID,
int  pixelID 
) const

Returns relative pixel efficiency (including CE, RQE and threshold efficiency)

Returns
pixel efficiency relative to nominal photocathode

Definition at line 223 of file TOPGeometryPar.cc.

224 {
225
226 auto channel = getChannelMapper().getChannel(pixelID);
227 auto pmtID = getChannelMapper().getPmtID(pixelID);
228
229 double RQE = getGeometry()->getPDETuningFactor(getPMTType(moduleID, pmtID));
230 if (m_channelRQE.isValid()) RQE *= m_channelRQE->getRQE(moduleID, channel);
231
232 double thrEffi = 1.0;
233 if (m_thresholdEff.isValid()) thrEffi = m_thresholdEff->getThrEff(moduleID, channel);
234
235 if (m_oldPayload) { // nominal QE is used
236 return RQE * thrEffi;
237 }
238
240
241 int id = getUniquePixelID(moduleID, pixelID);
242 return m_relEfficiencies[id] * RQE * thrEffi;
243 }
int getPmtID(int pixel) const
Returns PMT ID (1-based)
int getUniquePixelID(int moduleID, int pixelID) const
Returns unique pixel ID within the detector.
void prepareRelEfficiencies() const
Prepares a map of relative pixel efficiencies.
DBObjPtr< TOPCalChannelThresholdEff > m_thresholdEff
channel threshold effi.
double getPDETuningFactor(unsigned type) const
Returns photon detection efficiency tuning factor of a given PMT type.
Definition: TOPGeometry.cc:57

◆ getTTS()

const TOPNominalTTS & getTTS ( int  moduleID,
int  pmtID 
) const

Returns TTS of a PMT at given position.

Parameters
moduleIDslot ID
pmtIDPMT ID
Returns
TTS

Definition at line 255 of file TOPGeometryPar.cc.

256 {
257 auto pmtType = getPMTType(moduleID, pmtID);
258 return getGeometry()->getTTS(pmtType);
259 }
const TOPNominalTTS & getTTS(unsigned type) const
Returns time transition spread of a given PMT type.
Definition: TOPGeometry.cc:50

◆ getUniquePixelID()

int getUniquePixelID ( int  moduleID,
int  pixelID 
) const
inlineprivate

Returns unique pixel ID within the detector.

Parameters
moduleIDslot ID
pixelIDpixel ID
Returns
unique ID

Definition at line 270 of file TOPGeometryPar.h.

271 {
272 return (moduleID << 16) + pixelID;
273 }

◆ getUniquePmtID()

int getUniquePmtID ( int  moduleID,
int  pmtID 
) const
inlineprivate

Returns unique PMT ID within the detector.

Parameters
moduleIDslot ID
pmtIDPMT ID
Returns
unique ID

Definition at line 259 of file TOPGeometryPar.h.

260 {
261 return (moduleID << 16) + pmtID;
262 }

◆ Initialize() [1/2]

void Initialize ( )

Initialize from database.

Definition at line 89 of file TOPGeometryPar.cc.

90 {
91 m_fromDB = true;
92 m_valid = false;
93
94 if (m_geoDB) delete m_geoDB;
95 m_geoDB = new DBObjPtr<TOPGeometry>();
96
97 if (!m_geoDB->isValid()) {
98 B2ERROR("TOPGeometry: no payload found in database");
99 return;
100 }
101 if ((*m_geoDB)->getWavelengthFilter().getName().empty()) {
102 m_oldPayload = true;
103 B2WARNING("TOPGeometry: obsolete payload revision (pixel independent PDE) - please, check global tag");
104 }
105 if ((*m_geoDB)->getTTSes().empty()) {
106 B2WARNING("TOPGeometry: obsolete payload revision (nominal TTS only) - please, check global tag");
107 }
108 if ((*m_geoDB)->arePDETuningFactorsEmpty()) {
109 B2WARNING("TOPGeometry: old payload revision (before bugfix and update of optical properties)");
110 }
111
112 // Make sure that we abort as soon as the geometry changes
113 m_geoDB->addCallback([]() {
114 B2FATAL("Geometry cannot change during processing, "
115 "aborting (component TOP)");
116 });
117
119 if (!m_frontEndMapper.isValid()) {
120 B2ERROR("TOPFrontEndMaps: no payload found in database");
121 return;
122 }
123
126 B2ERROR("TOPChannelMaps: no payload found in database");
127 return;
128 }
129 m_valid = true;
130
132
133 }
void initialize(const GearDir &channelMapping)
Initialize from Gearbox (XML)
bool isValid() const
Checks if mapping is available.
Definition: ChannelMapper.h:80
void initialize(const GearDir &frontEndMapping)
Initialize from Gearbox (XML)
bool isValid() const
check if the mapping is available
void finalizeInitialization()
finalize initialization

◆ Initialize() [2/2]

void Initialize ( const GearDir content)

Initialize from Gearbox (XML)

Parameters
contentXML data directory

Definition at line 52 of file TOPGeometryPar.cc.

53 {
54
55 m_fromDB = false;
56 m_valid = false;
57
58 if (m_geo) delete m_geo;
59 m_geo = createConfiguration(content);
60 if (!m_geo->isConsistent()) {
61 B2ERROR("TOPGeometryPar::createConfiguration: geometry not consistently defined");
62 return;
63 }
64
65 GearDir frontEndMapping(content, "FrontEndMapping");
66 m_frontEndMapper.initialize(frontEndMapping);
68 return;
69 }
70
71 GearDir channelMapping0(content, "ChannelMapping[@type='IRS3B']");
72 m_channelMapperIRS3B.initialize(channelMapping0);
74 return;
75 }
76
77 GearDir channelMapping1(content, "ChannelMapping[@type='IRSX']");
78 m_channelMapperIRSX.initialize(channelMapping1);
80 return;
81 }
82 m_valid = true;
83
85
86 }
TOPGeometry * createConfiguration(const GearDir &content)
Create a parameter object from gearbox.
bool isConsistent() const override
Check for consistency of data members.
Definition: TOPGeometry.cc:126

◆ Instance()

TOPGeometryPar * Instance ( )
static

Static method to obtain the pointer to its instance.

Returns
pointer to the instance of this class.

Definition at line 43 of file TOPGeometryPar.cc.

44 {
45 if (!s_instance) {
47 }
48 return s_instance;
49 }
TOPGeometryPar()
Hidden constructor since it is a singleton class.

◆ integralOfQE()

double integralOfQE ( const std::vector< float > &  qe,
double  ce,
double  lambdaFirst,
double  lambdaStep 
) const
private

Returns integral of quantum efficiency over photon energies.

Parameters
qequantum efficiency data points
cecollection efficiency data points
lambdaFirstwavelenght of the first data point [nm]
lambdaStepwavelength step [nm]
Returns
integral [eV]

Definition at line 419 of file TOPGeometryPar.cc.

421 {
422 if (qe.empty()) return 0;
423
424 double s = 0;
425 double lambda = lambdaFirst;
426 double f1 = qe[0] / (lambda * lambda);
427 for (size_t i = 1; i < qe.size(); i++) {
428 lambda += lambdaStep;
429 double f2 = qe[i] / (lambda * lambda);
430 s += (f1 + f2) / 2;
431 f1 = f2;
432 }
433 return s * c_hc * lambdaStep * ce;
434 }

◆ isValid()

bool isValid ( ) const
inline

check if the geometry is available

Returns
true if available

Definition at line 65 of file TOPGeometryPar.h.

65{return m_valid;}

◆ mapPmtQEToPositions()

void mapPmtQEToPositions ( ) const
private

Maps PMT QE data to positions within the detector.

Definition at line 340 of file TOPGeometryPar.cc.

341 {
342 m_pmts.clear();
343
344 std::map<std::string, const TOPPmtQE*> map;
345 for (const auto& pmt : m_pmtQEData) {
346 map[pmt.getSerialNumber()] = &pmt;
347 }
348 for (const auto& pmt : m_pmtInstalled) {
349 int id = getUniquePmtID(pmt.getSlotNumber(), pmt.getPosition());
350 m_pmts[id] = map[pmt.getSerialNumber()];
351 }
352
353 B2INFO("TOPGeometryPar: QE of PMT's mapped to positions, size = " << m_pmts.size());
354 }

◆ mapPmtTypeToPositions()

void mapPmtTypeToPositions ( ) const
private

Maps PMT type to positions within the detector.

Definition at line 357 of file TOPGeometryPar.cc.

358 {
359 for (const auto& pmt : m_pmtInstalled) {
360 int id = getUniquePmtID(pmt.getSlotNumber(), pmt.getPosition());
361 m_pmtTypes[id] = pmt.getType();
362 }
363
364 B2INFO("TOPGeometryPar: PMT types mapped to positions, size = "
365 << m_pmtTypes.size());
366
367
368 std::set<unsigned> types;
369 for (const auto& pmt : m_pmtInstalled) {
370 types.insert(pmt.getType());
371 }
372 const auto* geo = getGeometry();
373 for (const auto& type : types) {
374 if (geo->getTTS(type).getPMTType() != type) {
375 B2WARNING("No TTS found for an installed PMT type. Nominal one will be used."
376 << LogVar("PMT type", type));
377 }
378 }
379
380 }

◆ prepareRelEfficiencies()

void prepareRelEfficiencies ( ) const
private

Prepares a map of relative pixel efficiencies.

Definition at line 383 of file TOPGeometryPar.cc.

384 {
385 m_relEfficiencies.clear();
386 if (m_pmts.empty()) mapPmtQEToPositions();
387
388 const auto* geo = getGeometry();
389
390 const auto& nominalQE = geo->getNominalQE();
391 double s0 = integralOfQE(nominalQE.getQE(), nominalQE.getCE(),
392 nominalQE.getLambdaFirst(), nominalQE.getLambdaStep());
393
394 for (const auto& module : geo->getModules()) {
395 auto moduleID = module.getModuleID();
396 const auto& pmtArray = module.getPMTArray();
397 int numPMTs = pmtArray.getSize();
398 int numPMTPixels = pmtArray.getPMT().getNumPixels();
399 for (int pmtID = 1; pmtID <= numPMTs; pmtID++) {
400 const auto* pmtQE = m_pmts[getUniquePmtID(moduleID, pmtID)];
401 for (int pmtPixel = 1; pmtPixel <= numPMTPixels; pmtPixel++) {
402 double s = 0;
403 if (pmtQE) {
404 s = integralOfQE(pmtQE->getQE(pmtPixel), pmtQE->getCE(m_BfieldOn),
405 pmtQE->getLambdaFirst(), pmtQE->getLambdaStep());
406 }
407 auto pixelID = pmtArray.getPixelID(pmtID, pmtPixel);
408 auto id = getUniquePixelID(moduleID, pixelID);
409 m_relEfficiencies[id] = s / s0;
410 }
411 }
412 }
413
414 B2INFO("TOPGeometryPar: pixel relative quantum efficiencies have been set, size = "
415 << m_relEfficiencies.size());
416 }
double integralOfQE(const std::vector< float > &qe, double ce, double lambdaFirst, double lambdaStep) const
Returns integral of quantum efficiency over photon energies.

◆ refractiveIndex()

double refractiveIndex ( double  lambda) const
private

Quartz refractive index (SellMeier equation)

Parameters
lambdaphoton wavelength [nm]
Returns
refractive index

Definition at line 959 of file TOPGeometryPar.cc.

960 {
961 // parameters of SellMeier equation (Matsuoka-san, 24.11.2018)
962 // from the specs of Corning HPFS 7980
963 // https://www.corning.com/media/worldwide/csm/documents/5bf092438c5546dfa9b08e423348317b.pdf
964 const double b[] = {0.683740494, 0.420323613, 0.585027480};
965 const double c[] = {0.00460352869, 0.0133968856, 64.4932732};
966
967 double x = pow(lambda * 0.001, 2);
968 double y = 1;
969 for (int i = 0; i < 3; i++) {
970 y += b[i] * x / (x - c[i]);
971 }
972 return sqrt(y);
973 }
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ setEnvelopeQE()

void setEnvelopeQE ( ) const
private

Constructs envelope of quantum efficiency from PMT data.

Definition at line 262 of file TOPGeometryPar.cc.

263 {
264 if (m_pmtQEData.getEntries() == 0) {
265 B2ERROR("DBArray TOPPmtQEs is empty");
266 return;
267 }
268
269 double lambdaFirst = 0;
270 for (const auto& pmt : m_pmtQEData) {
271 if (pmt.getLambdaFirst() > 0) {
272 lambdaFirst = pmt.getLambdaFirst();
273 break;
274 }
275 }
276 if (lambdaFirst == 0) {
277 B2ERROR("DBArray TOPPmtQEs: lambdaFirst of all PMT found to be less or equal 0");
278 return;
279 }
280 for (const auto& pmt : m_pmtQEData) {
281 if (pmt.getLambdaFirst() > 0) {
282 lambdaFirst = std::min(lambdaFirst, pmt.getLambdaFirst());
283 }
284 }
285
286 double lambdaStep = 0;
287 for (const auto& pmt : m_pmtQEData) {
288 if (pmt.getLambdaStep() > 0) {
289 lambdaStep = pmt.getLambdaStep();
290 break;
291 }
292 }
293 if (lambdaStep == 0) {
294 B2ERROR("DBArray TOPPmtQEs: lambdaStep of all PMT found to be less or equal 0");
295 return;
296 }
297 for (const auto& pmt : m_pmtQEData) {
298 if (pmt.getLambdaStep() > 0) {
299 lambdaStep = std::min(lambdaStep, pmt.getLambdaStep());
300 }
301 }
302
303 std::map<std::string, const TOPPmtInstallation*> map;
304 for (const auto& pmt : m_pmtInstalled) {
305 map[pmt.getSerialNumber()] = &pmt;
306 }
307 const auto* geo = getGeometry();
308
309 std::vector<float> envelopeQE;
310 for (const auto& pmt : m_pmtQEData) {
311 float ce = pmt.getCE(m_BfieldOn);
312 auto pmtInstalled = map[pmt.getSerialNumber()];
313 if (pmtInstalled) ce *= geo->getPDETuningFactor(pmtInstalled->getType());
314 if (pmt.getLambdaFirst() == lambdaFirst and pmt.getLambdaStep() == lambdaStep) {
315 const auto& envelope = pmt.getEnvelopeQE();
316 if (envelopeQE.size() < envelope.size()) {
317 envelopeQE.resize(envelope.size() - envelopeQE.size(), 0);
318 }
319 for (size_t i = 0; i < std::min(envelopeQE.size(), envelope.size()); i++) {
320 envelopeQE[i] = std::max(envelopeQE[i], envelope[i] * ce);
321 }
322 } else {
323 double lambdaLast = pmt.getLambdaLast();
324 int nExtra = (lambdaLast - lambdaFirst) / lambdaStep + 1 - envelopeQE.size();
325 if (nExtra > 0) envelopeQE.resize(nExtra, 0);
326 for (size_t i = 0; i < envelopeQE.size(); i++) {
327 float qe = pmt.getEnvelopeQE(lambdaFirst + lambdaStep * i);
328 envelopeQE[i] = std::max(envelopeQE[i], qe * ce);
329 }
330 }
331 }
332
333 m_envelopeQE.set(lambdaFirst, lambdaStep, 1.0, envelopeQE, "EnvelopeQE");
334
335 B2INFO("TOPGeometryPar: envelope of PMT dependent QE has been set");
336
337 }
void set(float lambdaFirst, float lambdaStep, float CE, const std::vector< float > &qe, const std::string &name)
Sets the object.
Definition: TOPNominalQE.h:55

Member Data Documentation

◆ c_hc

const double c_hc = 1239.84193
static

Planck constant times speed of light in [eV*nm].

Definition at line 175 of file TOPGeometryPar.h.

◆ m_BfieldOn

bool m_BfieldOn = true
private

true if B field is on

Definition at line 301 of file TOPGeometryPar.h.

◆ m_channelMapperIRS3B

ChannelMapper m_channelMapperIRS3B
private

channel-pixel mapper

Definition at line 306 of file TOPGeometryPar.h.

◆ m_channelMapperIRSX

ChannelMapper m_channelMapperIRSX
private

channel-pixel mapper

Definition at line 307 of file TOPGeometryPar.h.

◆ m_channelRQE

DBObjPtr<TOPCalChannelRQE> m_channelRQE
private

channel relative quantum effi.

Definition at line 313 of file TOPGeometryPar.h.

◆ m_envelopeQE

TOPNominalQE m_envelopeQE
mutableprivate

envelope quantum efficiency

Definition at line 317 of file TOPGeometryPar.h.

◆ m_fromDB

bool m_fromDB = false
private

parameters from database or Gearbox

Definition at line 298 of file TOPGeometryPar.h.

◆ m_frontEndMapper

FrontEndMapper m_frontEndMapper
private

front end electronics mapper

Definition at line 305 of file TOPGeometryPar.h.

◆ m_geo

TOPGeometry* m_geo = 0
private

geometry parameters from Gearbox

Definition at line 296 of file TOPGeometryPar.h.

◆ m_geoDB

DBObjPtr<TOPGeometry>* m_geoDB = 0
private

geometry parameters from database

Definition at line 297 of file TOPGeometryPar.h.

◆ m_oldPayload

bool m_oldPayload = false
private

true if old payload found in DB

Definition at line 300 of file TOPGeometryPar.h.

◆ m_pmtInstalled

OptionalDBArray<TOPPmtInstallation> m_pmtInstalled
private

PMT installation data.

Definition at line 311 of file TOPGeometryPar.h.

◆ m_pmtQEData

OptionalDBArray<TOPPmtQE> m_pmtQEData
private

quantum efficiencies

Definition at line 312 of file TOPGeometryPar.h.

◆ m_pmts

std::map<int, const TOPPmtQE*> m_pmts
mutableprivate

QE data mapped to positions.

Definition at line 318 of file TOPGeometryPar.h.

◆ m_pmtTypes

std::map<int, unsigned> m_pmtTypes
mutableprivate

PMT types mapped to positions.

Definition at line 320 of file TOPGeometryPar.h.

◆ m_relEfficiencies

std::map<int, double> m_relEfficiencies
mutableprivate

pixel relative QE

Definition at line 319 of file TOPGeometryPar.h.

◆ m_thresholdEff

DBObjPtr<TOPCalChannelThresholdEff> m_thresholdEff
private

channel threshold effi.

Definition at line 314 of file TOPGeometryPar.h.

◆ m_valid

bool m_valid = false
private

true if geometry is available

Definition at line 299 of file TOPGeometryPar.h.

◆ s_instance

TOPGeometryPar * s_instance = 0
staticprivate

Pointer to the class instance.

Definition at line 324 of file TOPGeometryPar.h.


The documentation for this class was generated from the following files: