File indexing completed on 2025-02-23 09:20:54
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
0027
0028
0029
0030
0031
0032
0033 #include "SteppingAction.hh"
0034
0035 #include "DetectorConstruction.hh"
0036 #include "HistoManager.hh"
0037 #include "RunAction.hh"
0038
0039 #include "G4ParticleTypes.hh"
0040 #include "G4RunManager.hh"
0041
0042 #include <G4RotationMatrix.hh>
0043 #include <G4ThreeVector.hh>
0044
0045
0046
0047 SteppingAction::SteppingAction(DetectorConstruction* det, RunAction* RuAct)
0048 : G4UserSteppingAction(), fDetector(det), fRunAction(RuAct)
0049 {}
0050
0051
0052
0053 SteppingAction::~SteppingAction() {}
0054
0055
0056
0057 void SteppingAction::UserSteppingAction(const G4Step* aStep)
0058 {
0059 G4StepPoint* prePoint = aStep->GetPreStepPoint();
0060
0061
0062 if (prePoint->GetTouchableHandle()->GetVolume() == fDetector->GetWorld()) return;
0063
0064
0065
0066
0067 G4RunManager::GetRunManager()->AbortEvent();
0068
0069
0070
0071 G4StepPoint* endPoint = aStep->GetPostStepPoint();
0072 G4String procName = endPoint->GetProcessDefinedStep()->GetProcessName();
0073 fRunAction->CountProcesses(procName);
0074
0075 if (procName == "msc" || procName == "muMsc" || procName == "stepMax") {
0076
0077
0078 G4ThreeVector position = endPoint->GetPosition();
0079 G4ThreeVector direction = endPoint->GetMomentumDirection();
0080
0081 G4double truePathLength = aStep->GetStepLength();
0082 G4double geomPathLength = position.x() + 0.5 * fDetector->GetBoxSize();
0083 G4double ratio = geomPathLength / truePathLength;
0084 fRunAction->SumPathLength(truePathLength, geomPathLength);
0085 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0086 analysisManager->FillH1(1, truePathLength);
0087 analysisManager->FillH1(2, geomPathLength);
0088 analysisManager->FillH1(3, ratio);
0089
0090 G4double yend = position.y(), zend = position.z();
0091 G4double lateralDisplacement = std::sqrt(yend * yend + zend * zend);
0092 fRunAction->SumLateralDisplacement(lateralDisplacement);
0093 analysisManager->FillH1(4, lateralDisplacement);
0094
0095 G4double psi = std::atan(lateralDisplacement / geomPathLength);
0096 fRunAction->SumPsi(psi);
0097 analysisManager->FillH1(5, psi);
0098
0099 G4double xdir = direction.x(), ydir = direction.y(), zdir = direction.z();
0100 G4double tetaPlane = std::atan2(ydir, xdir);
0101 fRunAction->SumTetaPlane(tetaPlane);
0102 analysisManager->FillH1(6, tetaPlane);
0103 tetaPlane = std::atan2(zdir, xdir);
0104 fRunAction->SumTetaPlane(tetaPlane);
0105 analysisManager->FillH1(6, tetaPlane);
0106
0107 G4double phiPos = std::atan2(zend, yend);
0108 analysisManager->FillH1(7, phiPos);
0109 G4double phiDir = std::atan2(zdir, ydir);
0110 analysisManager->FillH1(8, phiDir);
0111
0112 G4double phiCorrel = 0.;
0113 if (lateralDisplacement > 0.) phiCorrel = (yend * ydir + zend * zdir) / lateralDisplacement;
0114 fRunAction->SumPhiCorrel(phiCorrel);
0115 analysisManager->FillH1(9, phiCorrel);
0116 }
0117 else if (procName == "conv" || procName == "GammaToMuPair") {
0118
0119
0120 G4StepPoint* PrePoint = aStep->GetPreStepPoint();
0121 G4double EGamma = PrePoint->GetTotalEnergy();
0122 G4ThreeVector PGamma = PrePoint->GetMomentum();
0123 G4ThreeVector PolaGamma = PrePoint->GetPolarization();
0124
0125 G4double Eplus = -1;
0126 G4ThreeVector Pplus, Pminus, Precoil;
0127
0128 const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
0129
0130 const size_t Nsecondaries = (*secondary).size();
0131
0132
0133 if (Nsecondaries == 0) return;
0134
0135 for (size_t lp = 0; lp < std::min(Nsecondaries, size_t(2)); lp++) {
0136 if (((*secondary)[lp]->GetDefinition() == G4Electron::Definition())
0137 || ((*secondary)[lp]->GetDefinition() == G4MuonMinus::Definition()))
0138 {
0139 Pminus = (*secondary)[lp]->GetMomentum();
0140 }
0141 if (((*secondary)[lp]->GetDefinition() == G4Positron::Definition())
0142 || ((*secondary)[lp]->GetDefinition() == G4MuonPlus::Definition()))
0143 {
0144 Eplus = (*secondary)[lp]->GetTotalEnergy();
0145 Pplus = (*secondary)[lp]->GetMomentum();
0146 }
0147 }
0148
0149 if (Nsecondaries >= 3) {
0150 Precoil = (*secondary)[2]->GetMomentum();
0151 }
0152
0153 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0154
0155
0156
0157 G4ThreeVector z = PGamma.unit();
0158 G4ThreeVector x(1., 0., 0.);
0159
0160
0161
0162 if (PolaGamma.mag() != 0.0) {
0163 x = PolaGamma.unit();
0164 }
0165 else {
0166 x = z.orthogonal().unit();
0167 }
0168
0169 G4ThreeVector y = z;
0170 y = y.cross(x);
0171
0172 G4RotationMatrix GtoW(x, y, z);
0173 G4RotationMatrix WtoG = inverseOf(GtoW);
0174
0175 G4double angleE = Pplus.angle(Pminus) * EGamma;
0176 analysisManager->FillH1(10, angleE);
0177
0178 if (Nsecondaries >= 3) {
0179
0180 analysisManager->FillH1(11, std::log10(Precoil.mag()));
0181 analysisManager->FillH1(12, Precoil.transform(WtoG).phi());
0182 }
0183 G4double phiPlus = Pplus.transform(WtoG).phi();
0184 G4double phiMinus = Pminus.transform(WtoG).phi();
0185 analysisManager->FillH1(13, phiPlus);
0186 analysisManager->FillH1(14, std::cos(phiPlus + phiMinus) * -2.0);
0187 analysisManager->FillH1(15, Eplus / EGamma);
0188
0189 G4double phiPola = PolaGamma.transform(WtoG).phi();
0190 analysisManager->FillH1(16, phiPola);
0191 }
0192 }
0193
0194