Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAEventScheduler.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27#include <memory>
30#include "G4SystemOfUnits.hh"
31#include "G4UnitsTable.hh"
34#include "G4Timer.hh"
35#include "G4Scheduler.hh"
36#include "G4UserMeshAction.hh"
37#include "G4MoleculeCounter.hh"
39#include "G4Molecule.hh"
40
42 G4int pixel)
44 , fPixel(pixel)
45 , fInitialPixels(fPixel)
46 , fpMesh(new G4DNAMesh(boundingBox, fPixel))
47 , fpGillespieReaction(new G4DNAGillespieDirectMethod())
48 , fpEventSet(new G4DNAEventSet())
49 , fpUpdateSystem(new G4DNAUpdateSystemModel())
50{
51 if(!CheckingReactionRadius(fpMesh->GetResolution()))
52 {
53 G4String WarMessage = "resolution is not good : " +
54 std::to_string(fpMesh->GetResolution() / nm);
55 G4Exception("G4DNAEventScheduler::InitializeInMesh()", "WrongResolution",
56 JustWarning, WarMessage);
57 }
58}
59
61{
62 fCounterMap.clear();
63 if(fTimeToRecord.empty())
64 {
65 G4String WarMessage = "fTimeToRecord is empty ";
66 G4Exception("G4DNAEventScheduler::ClearAndReChargeCounter()",
67 "TimeToRecord is empty", JustWarning, WarMessage);
68 }
69 fLastRecoredTime = fTimeToRecord.begin();
70
71 if(G4VMoleculeCounter::Instance()->InUse()) // copy from MoleculeCounter
72 {
75 if(species.get() == nullptr)
76 {
77 return;
78 }
79 else if(species->empty())
80 {
82 return;
83 }
84 for(auto time_mol : fTimeToRecord)
85 {
86 if(time_mol > fStartTime)
87 {
88 continue;
89 }
90
91 for(auto molecule : *species)
92 {
94 molecule, time_mol);
95
96 if(n_mol < 0)
97 {
98 G4cerr << "G4DNAEventScheduler::ClearAndReChargeCounter() ::N "
99 "molecules not valid < 0 "
100 << G4endl;
101 G4Exception("", "N<0", FatalException, "");
102 }
103 fCounterMap[time_mol][molecule] = n_mol;
104 }
105 fLastRecoredTime++;
106 }
108 G4MoleculeCounter::Instance()->Use(false); // no more used
109 }
110 else
111 {
112 G4ExceptionDescription exceptionDescription;
113 exceptionDescription << "G4VMoleculeCounter is not used";
114 G4Exception("G4DNAEventScheduler::ClearAndReChargeCounter()",
115 "G4DNAEventScheduler010", JustWarning, exceptionDescription);
116 }
117}
118
119[[maybe_unused]] void G4DNAEventScheduler::AddTimeToRecord(const G4double& time)
120{
121 if(fTimeToRecord.find(time) == fTimeToRecord.end())
122 {
123 fTimeToRecord.insert(time);
124 }
125}
126
128
130{
131 auto pMainList = G4ITTrackHolder::Instance()->GetMainList();
132 std::map<G4VDNAMesh::Index, MapList> TrackKeyMap;
133 for(auto track : *pMainList)
134 {
135 auto molType = GetMolecule(track)->GetMolecularConfiguration();
136
137 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
139 if(pScavengerMaterial != nullptr &&
140 pScavengerMaterial->find(molType)) // avoid voxelize the scavenger
141 {
142 continue;
143 }
144
145 auto key = fpMesh->GetIndex(track->GetPosition());
146 if(TrackKeyMap.find(key) != TrackKeyMap.end())
147 {
148 std::map<MolType, size_t>& TrackTypeMap = TrackKeyMap[key];
149 if(TrackTypeMap.find(molType) != TrackTypeMap.end())
150 {
151 TrackTypeMap[molType]++;
152 }
153 else
154 {
155 TrackTypeMap[molType] = 1;
156 }
157 }
158 else
159 {
160 TrackKeyMap[key][molType] = 1;
161 }
162 }
163
164 for(auto& it : TrackKeyMap)
165 {
166 fpMesh->InitializeVoxel(it.first, std::move(it.second));
167 }
168}
169
171{
172 fPixel = pixel;
173 auto newMesh = new G4DNAMesh(fpMesh->GetBoundingBox(), fPixel);
174
175 auto begin = fpMesh->begin();
176 auto end = fpMesh->end();
177 std::map<G4VDNAMesh::Index, MapList> TrackKeyMap;
178 for(; begin != end; begin++)
179 {
180 auto index = std::get<0>(*begin);
181 auto newIndex = fpMesh->ConvertIndex(index, fPixel);
182 if(TrackKeyMap.find(newIndex) == TrackKeyMap.end())
183 {
184 TrackKeyMap[newIndex] = std::get<2>(*begin);
185 }
186 else
187 {
188 for(const auto& it : std::get<2>(*begin))
189 {
190 TrackKeyMap[newIndex][it.first] += it.second;
191 }
192 if(fVerbose > 1)
193 {
194 G4cout << " ReVoxelizing:: Old index : " << index
195 << " new index : " << fpMesh->ConvertIndex(index, fPixel)
196 << " number: " << std::get<2>(*begin).begin()->second << G4endl;
197 }
198 }
199 }
200 fpMesh.reset(newMesh);
201
202 for(auto& it : TrackKeyMap)
203 {
204 fpMesh->InitializeVoxel(it.first, std::move(it.second));
205 }
206}
208{
209 // find another solultion
210 fGlobalTime = fEndTime;
211
212 //
213 // RecordTime();//Last register for counter
214
215 if(fVerbose > 0)
216 {
217 G4cout << "End Processing and reset Gird, ScavengerTable, EventSet for new "
218 "simulation!!!!"
219 << G4endl;
220 }
221 fInitialized = false;
222 fTimeStep = 0;
223 fStepNumber = 0;
224 fGlobalTime = fStartTime;
225 fRunning = true;
226 fReactionNumber = 0;
227 fJumpingNumber = 0;
228
229 fpEventSet->RemoveEventSet();
230 fpMesh->Reset();
231}
232
234{
235 if(!fInitialized)
236 {
237 fPixel = fInitialPixels;
238 fpMesh = std::make_unique<G4DNAMesh>(fpMesh->GetBoundingBox(), fPixel);
239
240 // Scavenger();
241
242 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
244 if(pScavengerMaterial == nullptr)
245 {
246 G4cout << "There is no scavenger" << G4endl;
247 }
248 else
249 {
250 if(fVerbose > 1)
251 {
252 pScavengerMaterial->PrintInfo();
253 }
254 }
255
256 Voxelizing();
257 fpGillespieReaction->SetVoxelMesh(*fpMesh);
258 fpGillespieReaction->SetEventSet(fpEventSet.get());
259 fpGillespieReaction->SetTimeStep(
260 0); // reset fTimeStep = 0 in fpGillespieReaction
261 fpGillespieReaction->Initialize();
262 fpUpdateSystem->SetMesh(fpMesh.get());
264 fInitialized = true;
265 }
266
267 if(fVerbose > 0)
268 {
269 fpUpdateSystem->SetVerbose(1);
270 }
271
272 if(fVerbose > 2)
273 {
274 fpMesh->PrintMesh();
275 }
276}
278{
279 if(fPixel <= 1)
280 {
281 fRunning = false;
282 return;
283 }
284 // TEST /3
285 ReVoxelizing(fPixel / 2); //
286 // ReVoxelizing(fPixel/3);//
287
288 fpGillespieReaction->SetVoxelMesh(*fpMesh);
289 fpUpdateSystem->SetMesh(fpMesh.get());
290 fpGillespieReaction->Initialize();
291}
292
294{
295 if(fVerbose > 0)
296 {
297 G4cout
298 << "*** End Processing In Mesh and reset Mesh, EventSet for new Mesh!!!!"
299 << G4endl;
300 }
301 fpEventSet->RemoveEventSet();
302 fInitialized = false;
303 fIsChangeMesh = false;
304 fReactionNumber = 0;
305 fJumpingNumber = 0;
306 fStepNumberInMesh = 0;
307}
308
309G4double G4DNAEventScheduler::GetStartTime() const { return fStartTime; }
310
311G4double G4DNAEventScheduler::GetEndTime() const { return fEndTime; }
312
314{
315 return fTimeStep;
316}
317
318G4int G4DNAEventScheduler::GetVerbose() const { return fVerbose; }
319
321{
322 fMaxStep = max;
323}
324
326{
327 fStartTime = time;
328 fGlobalTime = fStartTime;
329}
330
331void G4DNAEventScheduler::Stop() { fRunning = false; }
333{
334 G4Timer localtimer;
335 if(fVerbose > 2)
336 {
337 localtimer.Start();
338 G4cout << "***G4DNAEventScheduler::Run*** for Pixel : " << fPixel << G4endl;
339 }
340 while(fEndTime > fGlobalTime && fRunning)
341 {
342 RunInMesh();
343 }
344 if(fVerbose > 2)
345 {
346 if(!fRunning)
347 {
348 G4cout << " StepNumber(" << fStepNumber << ") = MaxStep(" << fMaxStep
349 << ")" << G4endl;
350 }
351 else if(fEndTime <= fGlobalTime)
352 {
353 G4cout << " GlobalTime(" << fGlobalTime << ") > EndTime(" << fEndTime
354 << ")"
355 << " StepNumber : " << fStepNumber << G4endl;
356 }
357 localtimer.Stop();
358 G4cout << "***G4DNAEventScheduler::Ending::"
359 << G4BestUnit(fGlobalTime, "Time")
360 << " Events left : " << fpEventSet->size() << G4endl;
361 if(fVerbose > 1)
362 {
363 fpMesh->PrintMesh();
364 }
365 G4cout << " Computing Time : " << localtimer << G4endl;
366 }
367 Reset();
368}
369
371{
372 if(!fInitialized)
373 {
375 }
376 if(fVerbose > 0)
377 {
378 G4double resolution = fpMesh->GetResolution();
379 G4cout << "At Time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
380 << " the Mesh has " << fPixel << " x " << fPixel << " x " << fPixel
381 << " voxels with Resolution " << G4BestUnit(resolution, "Length")
382 << " during next "
383 << G4BestUnit(resolution * resolution * C / (6 * D), "Time")
384 << G4endl;
385 }
386
387 if(fVerbose > 2)
388 {
389 fpMesh->PrintMesh();
390 }
391
392 if(fpUserMeshAction != nullptr)
393 {
394 fpUserMeshAction->BeginOfMesh(fpMesh.get(), fGlobalTime);
395 }
396
397 // if diffusive jumping is avaiable, EventSet is never empty
398
399 while(!fpEventSet->Empty() && !fIsChangeMesh && fEndTime > fGlobalTime)
400 {
401 Stepping();
402 fGlobalTime = fTimeStep + fStartTime;
403
404 if(fpUserMeshAction != nullptr)
405 {
406 fpUserMeshAction->InMesh(fpMesh.get(), fGlobalTime);
407 }
408
409 if(fVerbose > 2)
410 {
411 G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime, "Time")
412 << " fTimeStep : " << G4BestUnit(fTimeStep, "Time") << G4endl;
413 }
414 G4double resolution = fpMesh->GetResolution();
415 fTransferTime = resolution * resolution * C / (6 * D);
416 if(fTransferTime == 0)
417 {
418 G4ExceptionDescription exceptionDescription;
419 exceptionDescription << "fTransferTime == 0";
420 G4Exception("G4DNAEventScheduler::RunInMesh", "G4DNAEventScheduler001",
421 FatalErrorInArgument, exceptionDescription);
422 }
423 if(fTransferTime < fTimeStep &&
424 fPixel != 1) // dont change Mesh if fPixel == 1
425 {
426 if(fSetChangeMesh)
427 {
428 if(fVerbose > 1)
429 {
430 G4cout << " Pixels : " << fPixel << " resolution : "
431 << G4BestUnit(fpMesh->GetResolution(), "Length")
432 << " fStepNumberInMesh : " << fStepNumberInMesh
433 << " at fGlobalTime : " << G4BestUnit(fGlobalTime, "Time")
434 << " at fTimeStep : " << G4BestUnit(fTimeStep, "Time")
435 << " fReactionNumber : " << fReactionNumber
436 << " transferTime : " << G4BestUnit(fTransferTime, "Time")
437 << G4endl;
438 }
439 fIsChangeMesh = true;
440 }
441 }
442 }
443
444 if(fVerbose > 1)
445 {
446 G4cout << "***G4DNAEventScheduler::Ending::"
447 << G4BestUnit(fGlobalTime, "Time")
448 << " Event left : " << fpEventSet->size() << G4endl;
449 G4cout << " Due to : ";
450 if(fpEventSet->Empty())
451 {
452 G4cout << "EventSet is Empty" << G4endl;
453 }
454 else if(fIsChangeMesh)
455 {
456 G4cout << "Changing Mesh from : " << fPixel
457 << " pixels to : " << fPixel / 2 << " pixels" << G4endl;
458 G4cout << "Info : ReactionNumber : " << fReactionNumber
459 << " JumpingNumber : " << fJumpingNumber << G4endl;
460 }
461 else if(fEndTime > fGlobalTime)
462 {
463 G4cout << " GlobalTime(" << fGlobalTime << ") > EndTime(" << fEndTime
464 << ")"
465 << " StepNumber : " << fStepNumber << G4endl;
466 }
467 if(fVerbose > 2)
468 {
469 fpMesh->PrintMesh();
470 }
471 G4cout << G4endl;
472 }
473
474 if(fpUserMeshAction != nullptr)
475 {
476 fpUserMeshAction->EndOfMesh(fpMesh.get(), fGlobalTime);
477 }
478 ResetInMesh();
479}
480
481void G4DNAEventScheduler::Stepping() // this event loop
482{
483 fStepNumber < fMaxStep ? fStepNumber++ : fRunning = false;
484 if(fpEventSet->size() > fpMesh->size())
485 {
486 G4ExceptionDescription exceptionDescription;
487 exceptionDescription
488 << "impossible that fpEventSet->size() > fpMesh->size()";
489 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler002",
490 FatalErrorInArgument, exceptionDescription);
491 }
492
493 auto selected = fpEventSet->begin();
494 auto index = (*selected)->GetIndex();
495
496 if(fVerbose > 1)
497 {
498 G4cout << "G4DNAEventScheduler::Stepping()*********************************"
499 "*******"
500 << G4endl;
501 (*selected)->PrintEvent();
502 }
503
504 // get selected time step
505 fTimeStep = (*selected)->GetTime();
506
507 // selected data
508 auto pJumping = (*selected)->GetJumpingData();
509 auto pReaction = (*selected)->GetReactionData();
510
511 fpUpdateSystem->SetGlobalTime(fTimeStep +
512 fStartTime); // this is just for printing
513 fpGillespieReaction->SetTimeStep(fTimeStep);
514 if(pJumping == nullptr && pReaction != nullptr)
515 {
516 fpUpdateSystem->UpdateSystem(index, *pReaction);
517 fpEventSet->RemoveEvent(selected);
518 // create new event
519 fpGillespieReaction->CreateEvent(index);
520 fReactionNumber++;
521 // recordTime in reaction
522 RecordTime();
523 }
524 else if(pJumping != nullptr && pReaction == nullptr)
525 {
526 // dont change this
527 fpUpdateSystem->UpdateSystem(index, *pJumping);
528 // save jumping Index before delete selected event
529 auto jumpingIndex = pJumping->second;
530 fpEventSet->RemoveEvent(selected);
531 // create new event
532 // should create Jumping before key
533 fpGillespieReaction->CreateEvent(jumpingIndex);
534 fpGillespieReaction->CreateEvent(index);
535 fJumpingNumber++;
536 }
537 else
538 {
539 G4ExceptionDescription exceptionDescription;
540 exceptionDescription << "pJumping == nullptr && pReaction == nullptr";
541 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler003",
542 FatalErrorInArgument, exceptionDescription);
543 }
544
545 if(fVerbose > 1)
546 {
547 G4cout << "G4DNAEventScheduler::Stepping::end "
548 "Print***********************************"
549 << G4endl;
550 G4cout << G4endl;
551 }
552 fStepNumberInMesh++;
553}
554
556{
557 fEndTime = endTime;
558}
559
561{
562 auto recordTime = *fLastRecoredTime;
563 if(fGlobalTime >= recordTime && fCounterMap[recordTime].empty())
564 {
565 auto begin = fpMesh->begin();
566 auto end = fpMesh->end();
567 for(; begin != end; begin++)
568 {
569 const auto& mapData = std::get<2>(*begin);
570 if(mapData.empty())
571 {
572 continue;
573 }
574 for(const auto& it : mapData)
575 {
576 fCounterMap[recordTime][it.first] += it.second;
577 }
578 }
579 fLastRecoredTime++;
580 }
581}
582
584{
585 G4cout << "CounterMap.size : " << fCounterMap.size() << G4endl;
586 for(const auto& i : fCounterMap)
587 {
588 auto map = i.second;
589 auto begin = map.begin(); //
590 auto end = map.end(); //
591 for(; begin != end; begin++)
592 {
593 auto molecule = begin->first;
594 auto number = begin->second;
595 if(number == 0)
596 {
597 continue;
598 }
599 G4cout << "molecule : " << molecule->GetName() << " number : " << number
600 << G4endl;
601 }
602 }
603}
604
607{
608 return fCounterMap;
609}
610
612 std::unique_ptr<G4UserMeshAction> pUserMeshAction)
613{
614 fpUserMeshAction = std::move(pUserMeshAction);
615}
616
617G4DNAMesh* G4DNAEventScheduler::GetMesh() const { return fpMesh.get(); }
618
619G4int G4DNAEventScheduler::GetPixels() const { return fPixel; }
620
622{
623 auto pMolecularReactionTable = G4DNAMolecularReactionTable::Instance();
624 auto reactionDataList = pMolecularReactionTable->GetVectorOfReactionData();
625 if(reactionDataList.empty())
626 {
627 G4cout << "reactionDataList.empty()" << G4endl;
628 return true;
629 }
630 else
631 {
632 for(auto it : reactionDataList)
633 {
634 if(it->GetEffectiveReactionRadius() >= resolution / CLHEP::pi)
635 {
636 G4cout << it->GetReactant1()->GetName() << " + "
637 << it->GetReactant2()->GetName() << G4endl;
638 G4cout << "G4DNAEventScheduler::ReactionRadius : "
639 << G4BestUnit(it->GetEffectiveReactionRadius(), "Length")
640 << G4endl;
641 G4cout << "resolution : " << G4BestUnit(resolution, "Length") << G4endl;
642 return false;
643 }
644 }
645 return true;
646 }
647}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4bool CheckingReactionRadius(G4double resolution)
G4DNAEventScheduler(const G4DNABoundingBox &boundingBox, G4int pixel)
void SetUserMeshAction(std::unique_ptr< G4UserMeshAction >)
std::map< G4double, MapCounter > GetCounterMap() const
void SetEndTime(const G4double &)
G4double GetTimeStep() const
G4double GetEndTime() const
~G4DNAEventScheduler() override
G4DNAMesh * GetMesh() const
G4double GetStartTime() const
void SetStartTime(G4double time)
std::map< MolType, G4int > MapCounter
void AddTimeToRecord(const G4double &time)
static G4DNAMolecularReactionTable * Instance()
G4TrackList * GetMainList(Key)
static G4ITTrackHolder * Instance()
std::unique_ptr< ReactantList > RecordedMolecules
static G4MoleculeCounter * Instance()
void ResetCounter() override
int GetNMoleculesAtTime(Reactant *molecule, double time)
RecordedMolecules GetRecordedMolecules()
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:530
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194
void Stop()
void Start()
void Use(G4bool flag=true)
static G4VMoleculeCounter * Instance()