12 #include "G4PrimaryParticle.hh" 13 #include "G4PrimaryVertex.hh" 14 #include "G4SDManager.hh" 51 label(iPSet.getUntrackedParameter<
std::
string>(
"instanceLabel",
"EcalValidInfo")) {
52 produces<PEcalValidInfo>(
label);
54 for (
int i = 0;
i < 26;
i++) {
76 for (
int i = 0;
i < 26;
i++) {
87 for (
int i = 0;
i < 26;
i++) {
185 for (
int i = 0;
i < 26;
i++) {
219 G4PrimaryParticle *thePrim =
nullptr;
220 int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
222 edm::LogWarning(
"EcalSimHitsValidProducer") <<
" No Vertex in this Event!";
224 for (
int i = 0;
i < nvertex;
i++) {
225 G4PrimaryVertex *avertex = (*evt)()->GetPrimaryVertex(
i);
226 if (avertex ==
nullptr)
227 edm::LogWarning(
"EcalSimHitsValidProducer") <<
" Pointer to vertex is NULL!";
229 float x0 = avertex->GetX0();
230 float y0 = avertex->GetY0();
231 float z0 = avertex->GetZ0();
232 float t0 = avertex->GetT0();
235 int npart = avertex->GetNumberOfParticle();
237 edm::LogWarning(
"EcalSimHitsValidProducer") <<
" No primary particle in this event";
239 if (thePrim ==
nullptr)
240 thePrim = avertex->GetPrimary(trackID);
247 if (thePrim !=
nullptr) {
248 double px = thePrim->GetPx();
249 double py = thePrim->GetPy();
250 double pz = thePrim->GetPz();
255 edm::LogWarning(
"EcalSimHitsValidProducer") <<
" Primary has p = 0 ; ";
260 thePID = thePrim->GetPDGcode();
262 edm::LogWarning(
"EcalSimHitsValidProducer") <<
" Could not find the primary particle!!";
266 G4HCofThisEvent *allHC = (*evt)()->GetHCofThisEvent();
267 int EBHCid = G4SDManager::GetSDMpointer()->GetCollectionID(
"EcalHitsEB");
268 int EEHCid = G4SDManager::GetSDMpointer()->GetCollectionID(
"EcalHitsEE");
269 int SEHCid = G4SDManager::GetSDMpointer()->GetCollectionID(
"EcalHitsES");
282 int theebhc_entries = theEBHC->entries();
283 for (
int j = 0;
j < theebhc_entries;
j++) {
290 float htheta = hpos.theta();
291 float heta = -
log(
tan(htheta * 0.5));
292 float hphi = hpos.phi();
305 MapType eemap, eezpmap, eezmmap;
306 int theeehc_entries = theEEHC->entries();
307 for (
int j = 0;
j < theeehc_entries;
j++) {
314 float htheta = hpos.theta();
315 float heta = -
log(
tan(htheta * 0.5));
316 float hphi = hpos.phi();
324 if (myEEid.
zside() == -1) {
329 if (myEEid.
zside() == 1) {
342 int thesehc_entries = theSEHC->entries();
343 for (
int j = 0;
j < thesehc_entries;
j++) {
348 if (esid.
zside() == -1) {
351 if (esid.
plane() == 1) {
354 }
else if (esid.
plane() == 2) {
359 if (esid.
zside() == 1) {
362 if (esid.
plane() == 1) {
365 }
else if (esid.
plane() == 2) {
374 if (eemap[eemaxid] > ebmap[ebmaxid]) {
377 int ix = myEEid.
ix();
378 int iy = myEEid.
iy();
379 int iz = myEEid.
zside();
395 int by = myEBid.
iphi();
396 int bz = myEBid.
zside();
412 G4StepPoint *preStepPoint = aStep->GetPreStepPoint();
413 const G4ThreeVector &hitPoint = preStepPoint->GetPosition();
414 G4VPhysicalVolume *currentPV = preStepPoint->GetPhysicalVolume();
415 const G4String &
name = currentPV->GetName();
417 crystal.assign(
name, 0, 4);
419 float Edeposit = aStep->GetTotalEnergyDeposit();
420 if (crystal ==
"EFRY" && Edeposit > 0.0) {
421 float z = hitPoint.z();
422 float detz = fabs(fabs(
z) - 3200);
423 int x0 = (
int)floor(detz / 8.9);
425 eEX0[x0] += Edeposit;
428 if (crystal ==
"EBRY" && Edeposit > 0.0) {
429 float x = hitPoint.x();
430 float y = hitPoint.y();
432 float detr =
r - 1290;
433 int x0 = (
int)floor(detr / 8.9);
435 eBX0[x0] += Edeposit;
441 int nCellInEta,
int nCellInPhi,
int CentralEta,
int CentralPhi,
int CentralZ,
MapType &fillmap,
MapType &themap) {
442 int goBackInEta = nCellInEta / 2;
443 int goBackInPhi = nCellInPhi / 2;
445 int startEta = CentralEta - goBackInEta;
446 int startPhi = CentralPhi - goBackInPhi;
449 for (
int ieta = startEta;
ieta < startEta + nCellInEta;
ieta++) {
450 for (
int iphi = startPhi;
iphi < startPhi + nCellInPhi;
iphi++) {
459 fillmap[
i++] = themap[
index];
464 if (fillmap[
i / 2] == themap[centerid])
471 int nCellInEta,
int nCellInPhi,
int CentralEta,
int CentralPhi,
int CentralZ,
MapType &fillmap,
MapType &themap) {
472 int goBackInEta = nCellInEta / 2;
473 int goBackInPhi = nCellInPhi / 2;
475 int startEta = CentralZ * CentralEta - goBackInEta;
476 int startPhi = CentralPhi - goBackInPhi;
479 for (
int ieta = startEta;
ieta < startEta + nCellInEta;
ieta++) {
480 for (
int iphi = startPhi;
iphi < startPhi + nCellInPhi;
iphi++) {
487 }
else if (
iphi > 360) {
492 fillmap[
i++] = themap[
index];
498 if (fillmap[
i / 2] == themap[ebcenterid])
506 float e012 = themap[0] + themap[1] + themap[2];
507 float e036 = themap[0] + themap[3] + themap[6];
508 float e678 = themap[6] + themap[7] + themap[8];
509 float e258 = themap[2] + themap[5] + themap[8];
511 if ((e012 > e678 || e012 == e678) && (e036 > e258 || e036 == e258))
512 E22 = themap[0] + themap[1] + themap[3] + themap[4];
513 else if ((e012 > e678 || e012 == e678) && (e036 < e258 || e036 == e258))
514 E22 = themap[1] + themap[2] + themap[4] + themap[5];
515 else if ((e012 < e678 || e012 == e678) && (e036 > e258 || e036 == e258))
516 E22 = themap[3] + themap[4] + themap[6] + themap[7];
517 else if ((e012 < e678 || e012 == e678) && (e036 < e258 || e036 == e258))
518 E22 = themap[4] + themap[5] + themap[7] + themap[8];
525 float e0_4 = themap[0] + themap[1] + themap[2] + themap[3] + themap[4];
526 float e0_20 = themap[0] + themap[5] + themap[10] + themap[15] + themap[20];
527 float e4_24 = themap[4] + themap[9] + themap[14] + themap[19] + themap[24];
528 float e0_24 = themap[20] + themap[21] + themap[22] + themap[23] + themap[24];
530 if ((e0_4 > e0_24 || e0_4 == e0_24) && (e0_20 > e4_24 || e0_20 == e4_24))
531 E44 = e33 + themap[0] + themap[1] + themap[2] + themap[3] + themap[5] + themap[10] + themap[15];
532 else if ((e0_4 > e0_24 || e0_4 == e0_24) && (e0_20 < e4_24 || e0_20 == e4_24))
533 E44 = e33 + themap[1] + themap[2] + themap[3] + themap[4] + themap[9] + themap[14] + themap[19];
534 else if ((e0_4 < e0_24 || e0_4 == e0_24) && (e0_20 > e4_24 || e0_20 == e4_24))
535 E44 = e33 + themap[5] + themap[10] + themap[15] + themap[20] + themap[21] + themap[22] + themap[23];
536 else if ((e0_4 < e0_24 || e0_4 == e0_24) && (e0_20 < e4_24 || e0_20 == e4_24))
537 E44 = e33 + themap[21] + themap[22] + themap[23] + themap[24] + themap[9] + themap[14] + themap[19];
543 int nCellInX,
int nCellInY,
int centralX,
int centralY,
int centralZ,
MapType &themap) {
545 float totalEnergy = 0.;
547 int goBackInX = nCellInX / 2;
548 int goBackInY = nCellInY / 2;
549 int startX = centralX - goBackInX;
550 int startY = centralY - goBackInY;
552 for (
int ix = startX;
ix < startX + nCellInX;
ix++) {
553 for (
int iy = startY;
iy < startY + nCellInY;
iy++) {
561 totalEnergy += themap[
index];
564 LogDebug(
"EcalSimHitsValidProducer")
565 <<
" EnergyInEEMatrix: ix - iy - E = " <<
ix <<
" " <<
iy <<
" " << themap[
index];
569 LogDebug(
"EcalSimHitsValidProducer") <<
" EnergyInEEMatrix: energy in " << nCellInX <<
" cells in x times " 570 << nCellInY <<
" cells in y matrix = " << totalEnergy <<
" for " << ncristals
577 int nCellInEta,
int nCellInPhi,
int centralEta,
int centralPhi,
int centralZ,
MapType &themap) {
579 float totalEnergy = 0.;
581 int goBackInEta = nCellInEta / 2;
582 int goBackInPhi = nCellInPhi / 2;
583 int startEta = centralZ * centralEta - goBackInEta;
584 int startPhi = centralPhi - goBackInPhi;
586 for (
int ieta = startEta;
ieta < startEta + nCellInEta;
ieta++) {
587 for (
int iphi = startPhi;
iphi < startPhi + nCellInPhi;
iphi++) {
594 }
else if (
iphi > 360) {
600 totalEnergy += themap[
index];
603 LogDebug(
"EcalSimHitsValidProducer")
604 <<
" EnergyInEBMatrix: ieta - iphi - E = " <<
ieta <<
" " <<
iphi <<
" " << themap[
index];
608 LogDebug(
"EcalSimHitsValidProducer") <<
" EnergyInEBMatrix: energy in " << nCellInEta <<
" cells in eta times " 609 << nCellInPhi <<
" cells in phi matrix = " << totalEnergy <<
" for " << ncristals
616 uint32_t unitWithMaxEnergy = 0;
619 MapType::iterator iter;
620 for (iter = themap.begin(); iter != themap.end(); iter++) {
623 unitWithMaxEnergy = (*iter).first;
627 LogDebug(
"EcalSimHitsValidProducer") <<
" Max energy of " <<
maxEnergy <<
" MeV was found in Unit id 0x" << std::hex
630 return unitWithMaxEnergy;
FloatVector eOfEEMinusCaloG4Hit
FloatVector eOfESCaloG4Hit
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
int ietaAbs() const
get the absolute value of the crystal ieta
FloatVector eOfEEPlusCaloG4Hit
FloatVector phiOfEECaloG4Hit
FloatVector eOfEEMinusCaloG4Hit
FloatVector tOfEECaloG4Hit
int iphi() const
get the crystal iphi
FloatVector etaOfEECaloG4Hit
FloatVector phiOfEBCaloG4Hit
math::XYZTLorentzVector theMomentum
FloatVector eOfESCaloG4Hit
bool fillEEMatrix(int nCellInEta, int nCellInPhi, int CentralEta, int CentralPhi, int CentralZ, MapType &fillmap, MapType &themap)
FloatVector etaOfEBCaloG4Hit
math::XYZTLorentzVector theVertex
math::XYZTLorentzVector theMomentum
FloatVector etaOfESCaloG4Hit
FloatVector eOfEBCaloG4Hit
math::XYZPoint getEntry() const
void fillEventInfo(PEcalValidInfo &)
uint32_t getUnitWithMaxEnergy(MapType &themap)
bool fillEBMatrix(int nCellInEta, int nCellInPhi, int CentralEta, int CentralPhi, int CentralZ, MapType &fillmap, MapType &themap)
FloatVector eOfEBCaloG4Hit
float energyInEBMatrix(int nCellInX, int nCellInY, int centralX, int centralY, int centralZ, MapType &themap)
Tan< T >::type tan(const T &t)
Abs< T >::type abs(const T &t)
FloatVector etaOfEBCaloG4Hit
FloatVector phiOfEECaloG4Hit
FloatVector phiOfESCaloG4Hit
FloatVector eOfEEPlusCaloG4Hit
std::map< uint32_t, float, std::less< uint32_t > > MapType
int zside() const
get the z-side of the crystal (1/-1)
FloatVector phiOfEBCaloG4Hit
float energyInEEMatrix(int nCellInX, int nCellInY, int centralX, int centralY, int centralZ, MapType &themap)
FloatVector tOfEECaloG4Hit
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
FloatVector tOfESCaloG4Hit
constexpr uint32_t rawId() const
get the raw id
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
XYZPointD XYZPoint
point in space with cartesian internal representation
double getEnergyDeposit() const
math::XYZTLorentzVector theVertex
uint32_t getUnitID() const
void produce(edm::Event &, const edm::EventSetup &) override
FloatVector etaOfESCaloG4Hit
FloatVector eOfEECaloG4Hit
float eCluster2x2(MapType &themap)
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
double getTimeSlice() const
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
~EcalSimHitsValidProducer() override
FloatVector tOfESCaloG4Hit
FloatVector phiOfESCaloG4Hit
Log< level::Warning, false > LogWarning
FloatVector etaOfEECaloG4Hit
EcalSimHitsValidProducer(const edm::ParameterSet &)
FloatVector tOfEBCaloG4Hit
FloatVector eOfEECaloG4Hit
FloatVector tOfEBCaloG4Hit
float eCluster4x4(float e33, MapType &themap)