File indexing completed on 2025-02-23 09:22:34
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026 #include "Par03Hit.hh"
0027
0028 #include "G4AttDef.hh"
0029 #include "G4AttDefStore.hh"
0030 #include "G4AttValue.hh"
0031 #include "G4Colour.hh"
0032 #include "G4LogicalVolume.hh"
0033 #include "G4SystemOfUnits.hh"
0034 #include "G4Tubs.hh"
0035 #include "G4UnitsTable.hh"
0036 #include "G4VVisManager.hh"
0037 #include "G4VisAttributes.hh"
0038
0039 G4ThreadLocal G4Allocator<Par03Hit>* Par03HitAllocator;
0040
0041 Par03Hit::Par03Hit() : G4VHit() {}
0042
0043
0044
0045 Par03Hit::~Par03Hit() = default;
0046
0047
0048
0049 Par03Hit::Par03Hit(const Par03Hit& aRight) : G4VHit()
0050 {
0051 fEdep = aRight.fEdep;
0052 fZId = aRight.fZId;
0053 fRhoId = aRight.fRhoId;
0054 fPhiId = aRight.fPhiId;
0055 fTime = aRight.fTime;
0056 fPos = aRight.fPos;
0057 fRot = aRight.fRot;
0058 fType = aRight.fType;
0059 fLogVol = aRight.fLogVol;
0060 }
0061
0062
0063
0064 const Par03Hit& Par03Hit::operator=(const Par03Hit& aRight)
0065 {
0066 fEdep = aRight.fEdep;
0067 fZId = aRight.fZId;
0068 fRhoId = aRight.fRhoId;
0069 fPhiId = aRight.fPhiId;
0070 fTime = aRight.fTime;
0071 fPos = aRight.fPos;
0072 fRot = aRight.fRot;
0073 fType = aRight.fType;
0074 fLogVol = aRight.fLogVol;
0075 return *this;
0076 }
0077
0078
0079
0080 int Par03Hit::operator==(const Par03Hit& aRight) const
0081 {
0082 return (fRhoId == aRight.fRhoId && fPhiId == aRight.fPhiId && fZId == aRight.fZId);
0083 }
0084
0085
0086
0087 void Par03Hit::Draw()
0088 {
0089 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
0090
0091 if (!pVVisManager->FilterHit(*this)) return;
0092
0093 if (fEdep < 0) return;
0094 if (pVVisManager) {
0095 G4Transform3D trans(fRot, fPos);
0096 G4VisAttributes attribs;
0097
0098 G4Tubs solid("draw", 0, 1 * cm, 1 * cm, 0, 0.05 * CLHEP::pi);
0099 if (fLogVol) {
0100 const G4VisAttributes* pVA = fLogVol->GetVisAttributes();
0101 if (pVA) attribs = *pVA;
0102
0103
0104 solid = *dynamic_cast<G4Tubs*>(fLogVol->GetSolid());
0105 double dR = solid.GetOuterRadius() - solid.GetInnerRadius();
0106 solid.SetInnerRadius(solid.GetInnerRadius() + fRhoId * dR);
0107 solid.SetOuterRadius(solid.GetOuterRadius() + fRhoId * dR);
0108 }
0109
0110 G4double colR = fType == 0 ? 0 : 1;
0111 G4double colG = fType == 0 ? 1 : 0;
0112 G4double colB = 0;
0113 G4Colour colour(colR, colG, colB, 0.5);
0114 attribs.SetColour(colour);
0115 attribs.SetForceSolid(true);
0116 pVVisManager->Draw(solid, attribs, trans);
0117 }
0118 }
0119
0120
0121
0122 const std::map<G4String, G4AttDef>* Par03Hit::GetAttDefs() const
0123 {
0124 G4bool isNew;
0125 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("Par03Hit", isNew);
0126 if (isNew) {
0127 (*store)["HitType"] = G4AttDef("HitType", "Hit Type", "Physics", "", "G4String");
0128 (*store)["Energy"] =
0129 G4AttDef("Energy", "Energy Deposited", "Physics", "G4BestUnit", "G4double");
0130 (*store)["Time"] = G4AttDef("Time", "Time", "Physics", "G4BestUnit", "G4double");
0131 (*store)["Pos"] = G4AttDef("Pos", "Position", "Physics", "G4BestUnit", "G4ThreeVector");
0132 }
0133 return store;
0134 }
0135
0136
0137
0138 std::vector<G4AttValue>* Par03Hit::CreateAttValues() const
0139 {
0140 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
0141 values->push_back(G4AttValue("HitType", "HadPar03Hit", ""));
0142 values->push_back(G4AttValue("Energy", G4BestUnit(fEdep, "Energy"), ""));
0143 values->push_back(G4AttValue("Time", G4BestUnit(fTime, "Time"), ""));
0144 values->push_back(G4AttValue("Pos", G4BestUnit(fPos, "Length"), ""));
0145 return values;
0146 }
0147
0148
0149
0150 void Par03Hit::Print()
0151 {
0152 std::cout << "\tHit " << fEdep / MeV << " MeV at " << fPos / cm << " cm (R,phi,z)= (" << fRhoId
0153 << ", " << fPhiId << ", " << fZId << "), " << fTime << " ns" << std::endl;
0154 }