File indexing completed on 2025-10-13 08:28:10
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 "Scorer.hh"
0027
0028 #include "InterPulseAction.hh"
0029 #include "PrimaryGeneratorAction.hh"
0030 #include "PulseAction.hh"
0031 #include "TimeStepAction.hh"
0032
0033 #include "G4AnalysisManager.hh"
0034 #include "G4DNAEventScheduler.hh"
0035 #include "G4DNAScavengerMaterial.hh"
0036 #include "G4Event.hh"
0037 #include "G4PhysicalConstants.hh"
0038 #include "G4RunManager.hh"
0039 #include "G4Scheduler.hh"
0040 #include "G4UIcmdWithADoubleAndUnit.hh"
0041 #include "G4UIcmdWithAnInteger.hh"
0042 #include "G4VChemistryWorld.hh"
0043
0044 #include <G4EventManager.hh>
0045 #include <G4MolecularConfiguration.hh>
0046 #include <G4SystemOfUnits.hh>
0047 #include <globals.hh>
0048
0049
0050
0051
0052 Dose::Dose()
0053 : G4UImessenger(),
0054 fpDoseDir(new G4UIdirectory("/scorer/Dose/")),
0055 fpAddDoseCutOff(new G4UIcmdWithADoubleAndUnit("/scorer/Dose/cutoff", this)),
0056 fpAddDoseToAbort(new G4UIcmdWithADoubleAndUnit("/scorer/Dose/abortedDose", this))
0057 {
0058 fpDoseDir->SetGuidance("Dose scorer commands");
0059 }
0060
0061
0062
0063 void Dose::SetNewValue(G4UIcommand* command, G4String newValue)
0064 {
0065 if (command == fpAddDoseCutOff.get()) {
0066 fDosesCutOff = G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(newValue);
0067 }
0068 if (command == fpAddDoseToAbort.get()) {
0069 fDosesToAbort = G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(newValue);
0070 }
0071 }
0072
0073
0074
0075 template<>
0076 void Scorer<Dose>::SetChemistryWorld(G4VChemistryWorld* pChemistryWorld)
0077 {
0078 fpChemistryWorld = pChemistryWorld;
0079 }
0080
0081
0082
0083 template<>
0084 G4VChemistryWorld* Scorer<Dose>::GetChemistryWorld() const
0085 {
0086 return fpChemistryWorld;
0087 }
0088
0089
0090
0091 template<>
0092 void Scorer<Dose>::clear()
0093 {
0094 fpScorer->fCumulatedDose = 0.;
0095 fpScorer->fPulseMax = 0;
0096 }
0097
0098 template<>
0099 void Scorer<Dose>::Initialize(G4HCofThisEvent* HCE)
0100 {
0101 clear();
0102 fpEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
0103 if (fHCID < 0) {
0104 fHCID = GetCollectionID(0);
0105 }
0106 HCE->AddHitsCollection(fHCID, (G4VHitsCollection*)fpEvtMap);
0107
0108 fPulseActionInfo =
0109 dynamic_cast<const InterPulseAction*>(G4RunManager::GetRunManager()->GetUserTrackingAction());
0110 }
0111
0112
0113
0114 template<>
0115 void Scorer<Dose>::EndOfEvent(G4HCofThisEvent*)
0116 {
0117 if (!G4RunManager::GetRunManager()->GetCurrentEvent()->IsAborted()) {
0118 fpEvtMap->add(0, fpScorer->fDosesCutOff);
0119 }
0120 fpScorer->fCumulatedDose = 0.;
0121 fpScorer->fPulseMax = 0;
0122 }
0123
0124
0125
0126 template<>
0127 G4bool Scorer<Dose>::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0128 {
0129 auto currentEvent = G4EventManager::GetEventManager();
0130 const G4Track* track = aStep->GetTrack();
0131 auto boundingBox = fpChemistryWorld->GetChemistryBoundary();
0132
0133 G4double V = boundingBox->Volume() / cm3;
0134 G4double edep = aStep->GetTotalEnergyDeposit();
0135 if (edep == 0.) {
0136 return false;
0137 }
0138 G4double eToGray = (0.001 * V * 6.242e+18);
0139 if (eToGray != 0) {
0140 (fpScorer->fCumulatedDose) += ((edep) / eV) / eToGray;
0141 }
0142 if (track->GetParentID() == 0 && aStep->IsFirstStepInVolume()) {
0143 auto pulseInfo = dynamic_cast<PulseInfo*>(track->GetUserInformation());
0144 if (pulseInfo != nullptr) {
0145 auto delayedTime = pulseInfo->GetDelayedTime();
0146 if (fpScorer->fPulseMax < delayedTime) {
0147 fpScorer->fPulseMax = delayedTime;
0148 }
0149 }
0150 G4double DoseInGray = fpScorer->fCumulatedDose;
0151 if (DoseInGray > fpScorer->fDosesCutOff / gray) {
0152 G4cout << "_____________________________________________________________________________"
0153 << G4endl;
0154 auto name = currentEvent->GetConstCurrentEvent()
0155 ->GetPrimaryVertex()
0156 ->GetPrimary()
0157 ->GetParticleDefinition()
0158 ->GetParticleName();
0159 auto energy =
0160 currentEvent->GetConstCurrentEvent()->GetPrimaryVertex()->GetPrimary()->GetKineticEnergy();
0161
0162 G4cout << "Beam line : " << "(" << name << ", " << energy << " MeV)" << G4endl;
0163 G4cout << "Cut-off dose : " << fpScorer->fDosesCutOff / gray << " Gy" << G4endl;
0164 G4cout << "Stop at actual dose : " << DoseInGray << " Gy" << G4endl;
0165 if (fPulseActionInfo != nullptr) {
0166 auto numberOfPulse = fPulseActionInfo->GetNumberOfPulse();
0167 auto DIT = fPulseActionInfo->GetPulsePeriod();
0168 G4cout << "DIT : " << DIT / CLHEP::ms << " ms" << G4endl;
0169 G4cout << "Pulse number : " << numberOfPulse << G4endl;
0170 }
0171
0172 const auto generatorAction = static_cast<const PrimaryGeneratorAction*>(
0173 G4RunManager::GetRunManager()->GetUserPrimaryGeneratorAction());
0174 const auto particleGun = generatorAction->GetSPGun();
0175 auto NumberOfParticlesGeneratedinOneEvent = particleGun->GetNumberOfParticles();
0176
0177 G4cout << "Beam duration : " << fpScorer->fPulseMax / second << " s" << G4endl;
0178 if(fpScorer->fPulseMax != 0){
0179 G4cout << "Actual dose rate : " << DoseInGray / (fpScorer->fPulseMax / second) << " Gy/s"
0180 << G4endl;
0181 }else{
0182 G4cout << "Actual dose rate : " << "infinite"
0183 << G4endl;
0184 }
0185
0186 G4cout << "Track number : "
0187 << NumberOfParticlesGeneratedinOneEvent
0188 - currentEvent->GetStackManager()->GetNUrgentTrack() - 1
0189 << " tracks" << G4endl;
0190 G4cout << "Irradiated volume : " << V * 1e+12
0191 << " um3 (" << 2 * boundingBox->halfSideLengthInX() / um << " x "
0192 << 2 * boundingBox->halfSideLengthInY() / um << " x "
0193 << 2 * boundingBox->halfSideLengthInZ() / um << ")" << G4endl;
0194 auto ChemComponent = fpChemistryWorld->GetChemicalComponent();
0195 for (const auto& it : ChemComponent) {
0196 if (fH2O == it.first) continue;
0197 if (fOHm == it.first) continue;
0198 if (it.first == fH3Op) {
0199 G4cout << "pH : " << -std::log10(it.second * (mole * liter)) << G4endl;
0200 continue;
0201 }
0202 G4cout << it.first->GetName() << " : " << it.second * (mole * liter) << " M "
0203 << G4endl;
0204 }
0205
0206 G4cout << "Total deposit energy : " << fpScorer->fCumulatedDose * eToGray << " eV" << G4endl;
0207 G4double DoseAbort;
0208 if (fpScorer->fDosesToAbort == 0) {
0209 DoseAbort = fpScorer->fDosesCutOff / gray + 0.5 * fpScorer->fDosesCutOff / gray;
0210 }
0211 else {
0212 DoseAbort = fpScorer->fDosesToAbort / gray;
0213 }
0214
0215 G4cout << "Dose to abort : " << DoseAbort << " Gy" << G4endl;
0216 if (DoseInGray > DoseAbort) {
0217 G4cout << "Abort this beam line (" << name << ", " << energy
0218 << " MeV) at actual dose: " << DoseInGray << " Gy" << G4endl;
0219 G4RunManager::GetRunManager()->AbortEvent();
0220 }
0221 G4cout << "_____________________________________________________________________________"
0222 << G4endl;
0223 auto myTrack = ((G4Track*)track);
0224 myTrack->SetTrackStatus(fStopAndKill);
0225 auto secondaries = track->GetStep()->GetSecondaryInCurrentStep();
0226 if (!secondaries->empty()) {
0227 for (auto it : *(secondaries)) {
0228 if (it != nullptr) {
0229 ((G4Track*)it)->SetTrackStatus(fStopAndKill);
0230 }
0231 }
0232 }
0233 currentEvent->GetStackManager()->ClearUrgentStack();
0234 }
0235 }
0236 return true;
0237 }
0238
0239
0240
0241
0242 Gvalues::Gvalues()
0243 : G4UImessenger(),
0244 fTimeLimit(G4Scheduler::Instance()->GetEndTime()),
0245 fSpeciesdir(new G4UIdirectory("/scorer/Gvalues/")),
0246 fTimeBincmd(new G4UIcmdWithAnInteger("/scorer/Gvalues/nOfTimeBins", this)),
0247 fAddTimeToRecordcmd(new G4UIcmdWithADoubleAndUnit("/scorer/Gvalues/addTimeToRecord", this))
0248 {
0249 fSpeciesdir->SetGuidance("ScoreSpecies commands");
0250 }
0251
0252
0253
0254 void Gvalues::SetNewValue(G4UIcommand* command, G4String newValue)
0255 {
0256 if (command == fAddTimeToRecordcmd.get()) {
0257 G4double cmdTime = fAddTimeToRecordcmd->GetNewDoubleValue(newValue);
0258 if (fTimeLimit >= cmdTime) {
0259 AddTimeToRecord(cmdTime);
0260 }
0261 else {
0262 AddTimeToRecord(fTimeLimit);
0263 }
0264 }
0265 if (command == fTimeBincmd.get()) {
0266 G4int cmdBins = fTimeBincmd->GetNewIntValue(newValue);
0267 G4double timeMin = 1 * ps;
0268 G4double timeMax = G4Scheduler::Instance()->GetEndTime() - timeMin;
0269 G4double timeLogMin = std::log10(timeMin);
0270 G4double timeLogMax = std::log10(timeMax);
0271 for (int i = 0; i <= cmdBins; i++) {
0272 AddTimeToRecord(std::pow(10, timeLogMin + i * (timeLogMax - timeLogMin) / cmdBins));
0273 }
0274 }
0275 }
0276
0277
0278 void Gvalues::WriteInfo(G4VAnalysisManager* analysisManager, const std::string& out)
0279 {
0280 G4int NtupleID = analysisManager->CreateNtuple("info", "Simulation");
0281 analysisManager->CreateNtupleDColumn(NtupleID, "Dose");
0282 analysisManager->FinishNtuple(NtupleID);
0283 analysisManager->FillNtupleDColumn(NtupleID, 0, std::stod(out));
0284 analysisManager->AddNtupleRow(NtupleID);
0285 }
0286
0287 void Gvalues::WriteGvalues(G4VAnalysisManager* analysisManager)
0288 {
0289 G4int NtupleID = analysisManager->CreateNtuple("Gvalue", "Gvalue");
0290 analysisManager->CreateNtupleIColumn(NtupleID, "speciesID");
0291 analysisManager->CreateNtupleIColumn(NtupleID, "number");
0292 analysisManager->CreateNtupleIColumn(NtupleID, "nEvent");
0293 analysisManager->CreateNtupleSColumn(NtupleID, "speciesName");
0294 analysisManager->CreateNtupleDColumn(NtupleID, "time");
0295 analysisManager->CreateNtupleDColumn(NtupleID, "sumG");
0296 analysisManager->CreateNtupleDColumn(NtupleID, "sumG2");
0297 analysisManager->CreateNtupleDColumn(NtupleID, "TotalDose");
0298 analysisManager->CreateNtupleDColumn(NtupleID, "TotalDose2");
0299 analysisManager->CreateNtupleDColumn(NtupleID, "TotalDoseRate");
0300 analysisManager->CreateNtupleDColumn(NtupleID, "TotalDoseRate2");
0301 analysisManager->FinishNtuple(NtupleID);
0302
0303 for (const auto& it_map1 : fSpeciesInfoPerTime) {
0304 const InnerSpeciesMap& map2 = it_map1.second;
0305 for (auto it_map2 : map2) {
0306 double time = it_map1.first;
0307 auto species = it_map2.first;
0308 const G4String& name = species->GetName();
0309 int molID = it_map2.first->GetMoleculeID();
0310 auto number = it_map2.second.fNumber;
0311 double G = it_map2.second.fG;
0312 double G2 = it_map2.second.fG2;
0313
0314 analysisManager->FillNtupleIColumn(NtupleID, 0, molID);
0315 analysisManager->FillNtupleIColumn(NtupleID, 1, number);
0316 analysisManager->FillNtupleIColumn(NtupleID, 2, fNEvent);
0317 analysisManager->FillNtupleSColumn(NtupleID, 3, name);
0318 analysisManager->FillNtupleDColumn(NtupleID, 4, time);
0319 analysisManager->FillNtupleDColumn(NtupleID, 5, G);
0320 analysisManager->FillNtupleDColumn(NtupleID, 6, G2);
0321 analysisManager->FillNtupleDColumn(NtupleID, 7, fTotalDose);
0322 analysisManager->FillNtupleDColumn(NtupleID, 8, fTotalDose2);
0323 analysisManager->FillNtupleDColumn(NtupleID, 9, fTotalDoseRate);
0324 analysisManager->FillNtupleDColumn(NtupleID, 10, fTotalDoseRate2);
0325 analysisManager->AddNtupleRow(NtupleID);
0326 }
0327 }
0328 }
0329
0330
0331 void Gvalues::WriteWithAnalysisManager(G4VAnalysisManager* analysisManager, const std::string& out)
0332 {
0333 WriteInfo(analysisManager, out);
0334 WriteGvalues(analysisManager);
0335
0336 }
0337
0338
0339 template<>
0340 void Scorer<Gvalues>::clear()
0341 {
0342 fpEvtMap->clear();
0343 fpScorer->fNEvent = 0;
0344 fpScorer->fEdep = 0;
0345 fpScorer->fSpeciesInfoPerTime.clear();
0346 }
0347
0348
0349
0350 template<>
0351 void Scorer<Gvalues>::Initialize(G4HCofThisEvent* HCE)
0352 {
0353 fpEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
0354 if (fHCID < 0) {
0355 fHCID = GetCollectionID(0);
0356 }
0357 HCE->AddHitsCollection(fHCID, (G4VHitsCollection*)fpEvtMap);
0358 }
0359
0360
0361
0362 template<>
0363 void Scorer<Gvalues>::SetChemistryWorld(G4VChemistryWorld* pChemistryWorld)
0364 {
0365 fpChemistryWorld = pChemistryWorld;
0366 }
0367
0368
0369
0370 template<>
0371 G4VChemistryWorld* Scorer<Gvalues>::GetChemistryWorld() const
0372 {
0373 return fpChemistryWorld;
0374 }
0375
0376
0377
0378 template<>
0379 G4bool Scorer<Gvalues>::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0380 {
0381 const G4Track* track = aStep->GetTrack();
0382 auto pulseInfo = dynamic_cast<PulseInfo*>(track->GetUserInformation());
0383 if (pulseInfo != nullptr) {
0384 auto delayedTime = pulseInfo->GetDelayedTime();
0385 if (fpScorer->fPulseMax < delayedTime) {
0386 fpScorer->fPulseMax = delayedTime;
0387 }
0388 }
0389
0390 G4double edep = aStep->GetTotalEnergyDeposit();
0391 if (edep == 0.) {
0392 return FALSE;
0393 }
0394 edep *= aStep->GetPreStepPoint()->GetWeight();
0395 G4int index = GetIndex(aStep);
0396 fpEvtMap->add(index, edep);
0397 (fpScorer->fEdep) += edep;
0398 return TRUE;
0399 }
0400
0401
0402
0403 template<>
0404 void Scorer<Gvalues>::SaveScavengerChange()
0405 {
0406 auto pScavengerMaterial =
0407 dynamic_cast<G4DNAScavengerMaterial*>(G4Scheduler::Instance()->GetScavengerMaterial());
0408 if (pScavengerMaterial == nullptr) {
0409 G4ExceptionDescription errMsg;
0410 errMsg << "pScavengerMaterial == nullptr";
0411 G4Exception("Scorer<Gvalues>::SaveScavengerChange()", "SaveScavengerChange",
0412 FatalErrorInArgument, errMsg);
0413 }
0414 auto scavengerList = pScavengerMaterial->GetScavengerList();
0415 auto V = fpChemistryWorld->GetChemistryBoundary()->Volume();
0416 auto iniC = fpChemistryWorld->GetChemicalComponent();
0417
0418 G4double eToGray = (0.001 * V * 6.242e+18);
0419 if (eToGray != 0) {
0420 G4double dose = 1000 * (fpScorer->fEdep / eV) / eToGray;
0421 for (const auto& it : scavengerList) {
0422 if (it == fH2O || fH3Op == it || fOHm == it) {
0423 continue;
0424 }
0425 for (auto time_mol : fpScorer->fTimeToRecord) {
0426 int64_t n_mol = pScavengerMaterial->GetNMoleculesAtTime(it, time_mol);
0427 if (n_mol < 0) {
0428 G4ExceptionDescription errMsg;
0429 errMsg << "SaveScavengerChange()::N molecules not valid < 0 : " << it->GetName()
0430 << " N : " << n_mol << G4endl;
0431 G4Exception("", "N<0", FatalException, errMsg);
0432 }
0433
0434 Gvalues::SpeciesInfo& molInfo = fpScorer->fSpeciesInfoPerTime[time_mol][it];
0435 molInfo.fNumber += n_mol;
0436 if (V > 0) {
0437 auto concentration = n_mol / (Avogadro * V );
0438
0439
0440 auto inuM = (iniC[it] - concentration) * (mole * liter) * 1e6;
0441
0442 if (dose > 0) {
0443 G4double gValue = inuM / dose;
0444 molInfo.fG += gValue;
0445 molInfo.fG2 += gValue * gValue;
0446 }
0447 }
0448 }
0449 }
0450 }
0451 }
0452
0453
0454
0455 template<>
0456 void Scorer<Gvalues>::SaveMoleculeCounter()
0457 {
0458 if (fpEventScheduler == nullptr) {
0459 G4Exception("fpEventScheduler == nullptr", "Scorer<Gvalues>::SaveMoleculeCounter()",
0460 FatalException, "fpEventScheduler == nullptr");
0461 }
0462 else {
0463 auto counterMap = fpEventScheduler->GetCounterMap();
0464 if (counterMap.empty()) {
0465 G4Exception("No counter", "Scorer<Gvalues>::SaveMoleculeCounter()", JustWarning,
0466 "CounterMap is not used");
0467 return;
0468 }
0469 for (const auto& map_mol : counterMap) {
0470 auto time_mol = map_mol.first;
0471 for (auto it_mol : map_mol.second) {
0472 auto molecule = it_mol.first;
0473 if (molecule == fO2) {
0474 continue;
0475 }
0476 int n_mol = it_mol.second;
0477
0478 if (n_mol < 0) {
0479 G4ExceptionDescription errMsg;
0480 errMsg << "N molecules not valid < 0 " << " molecule : " << it_mol.first->GetName()
0481 << " N : " << n_mol << G4endl;
0482 G4Exception("", "N<0", FatalException, errMsg);
0483 }
0484
0485 if (fpScorer->fEdep > 0) {
0486 Gvalues::SpeciesInfo& molInfo = fpScorer->fSpeciesInfoPerTime[time_mol][molecule];
0487 molInfo.fNumber += n_mol;
0488 G4double gValue = (n_mol / (fpScorer->fEdep / eV)) * 100.;
0489
0490 molInfo.fG += gValue;
0491 molInfo.fG2 += gValue * gValue;
0492 }
0493 }
0494 }
0495
0496 auto boundingBox = fpChemistryWorld->GetChemistryBoundary();
0497 G4double V = boundingBox->Volume() / cm3;
0498 G4double eToGray = (0.001 * V * 6.242e+18);
0499 if (eToGray != 0) {
0500 G4double dose = (fpScorer->fEdep / eV) / eToGray;
0501 if (fpScorer->fPulseMax != 0) {
0502 G4double doseRate = dose / (fpScorer->fPulseMax / second);
0503 fpScorer->fTotalDose += dose;
0504 fpScorer->fTotalDose2 += dose * dose;
0505 fpScorer->fTotalDoseRate += doseRate;
0506 fpScorer->fTotalDoseRate2 += doseRate * doseRate;
0507 }
0508 }
0509 }
0510 }
0511
0512
0513
0514 template<>
0515 void Scorer<Gvalues>::EndOfEvent(G4HCofThisEvent*)
0516 {
0517 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted()) {
0518 fpScorer->fEdep = 0.;
0519 fpEventScheduler->ResetCounter();
0520 return;
0521 }
0522
0523 SaveScavengerChange();
0524 SaveMoleculeCounter();
0525
0526 ++(fpScorer->fNEvent);
0527 fpScorer->fEdep = 0.;
0528
0529 fpEventScheduler->ResetCounter();
0530 }
0531
0532
0533
0534 template<>
0535 void Scorer<Gvalues>::AbsorbResultsFromWorkerScorer(G4VPrimitiveScorer* workerScorer)
0536 {
0537 auto right = dynamic_cast<Scorer<Gvalues>*>(dynamic_cast<G4VPrimitiveScorer*>(workerScorer));
0538
0539 if (right == nullptr) {
0540 return;
0541 }
0542 if (right == this) {
0543 return;
0544 }
0545 auto it_map1 = right->fpScorer->fSpeciesInfoPerTime.begin();
0546 auto end_map1 = right->fpScorer->fSpeciesInfoPerTime.end();
0547
0548 for (; it_map1 != end_map1; ++it_map1) {
0549 Gvalues::InnerSpeciesMap& map2 = it_map1->second;
0550 auto it_map2 = map2.begin();
0551 auto end_map2 = map2.end();
0552
0553 for (; it_map2 != end_map2; ++it_map2) {
0554 Gvalues::SpeciesInfo& molInfo = fpScorer->fSpeciesInfoPerTime[it_map1->first][it_map2->first];
0555 molInfo.fNumber += it_map2->second.fNumber;
0556 molInfo.fG += it_map2->second.fG;
0557 molInfo.fG2 += it_map2->second.fG2;
0558 }
0559 }
0560
0561
0562 fpScorer->fTotalDose += right->fpScorer->fTotalDose;
0563 fpScorer->fTotalDose2 += right->fpScorer->fTotalDose2;
0564 right->fpScorer->fTotalDose = 0;
0565 right->fpScorer->fTotalDose2 = 0.;
0566
0567
0568 fpScorer->fTotalDoseRate += right->fpScorer->fTotalDoseRate;
0569 fpScorer->fTotalDoseRate2 += right->fpScorer->fTotalDoseRate2;
0570 right->fpScorer->fTotalDoseRate = 0;
0571 right->fpScorer->fTotalDoseRate2 = 0.;
0572
0573 right->fpScorer->fSpeciesInfoPerTime.clear();
0574 fpScorer->fNEvent += right->fpScorer->fNEvent;
0575 right->fpScorer->fNEvent = 0;
0576 right->fpScorer->fEdep = 0.;
0577 }
0578
0579
0580
0581 template<>
0582 void Scorer<Gvalues>::OutputAndClear(const std::string& dose)
0583 {
0584 if (G4Threading::IsWorkerThread()) {
0585 return;
0586 }
0587 G4VAnalysisManager* analysisManager = G4AnalysisManager::Instance();
0588 if (analysisManager != nullptr) {
0589 this->fpScorer->WriteWithAnalysisManager(analysisManager, dose);
0590 }
0591 fpScorer->fNEvent = 0;
0592 fpScorer->fSpeciesInfoPerTime.clear();
0593 }
0594
0595