Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/utils/include/G4DNAEventSet.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 #ifndef G4DNAEventSet_hh
 28 #define G4DNAEventSet_hh 1
 29 #include <list>
 30 #include <map>
 31 #include <G4memory.hh>
 32 #include <set>
 33 #include "G4DNAMesh.hh"
 34 #include <variant>
 35 class G4DNAMolecularReactionTable;
 36 class G4DNAMolecularReactionData;
 37 class Event
 38 {
 39  public:
 40   using Index        = G4VDNAMesh::Index;
 41   using MolType      = const G4MolecularConfiguration*;
 42   using JumpingData  = std::pair<MolType, Index>;
 43   using ReactionData = const G4DNAMolecularReactionData;
 44 
 45   // to test C++17
 46   // using Data = std::variant<std::unique_ptr<JumpingData>, ReactionData*>
 47   using Data = std::pair<std::unique_ptr<JumpingData>, ReactionData*>;
 48 
 49   Event(const G4double& time, const Index& index, ReactionData*);
 50   Event(const G4double& time, const Index& index, std::unique_ptr<JumpingData>&&);
 51 
 52   virtual ~Event();
 53   inline G4double GetTime() const { return fTimeStep; }
 54   inline Index GetIndex() const { return fIndex; }
 55   void PrintEvent() const;
 56   JumpingData* GetJumpingData() const { return std::get<0>(fData).get(); }
 57   ReactionData* GetReactionData() const { return std::get<1>(fData); }
 58 
 59  private:
 60   G4double fTimeStep;
 61   Index fIndex;
 62   Data fData;
 63 };
 64 
 65 struct comparatorEventSet
 66 {
 67   G4bool operator()(std::unique_ptr<Event> const& rhs,
 68                     std::unique_ptr<Event> const& lhs) const;
 69 };
 70 
 71 class IEventSet
 72 {
 73  public:
 74   IEventSet()  = default;
 75   ~IEventSet() = default;
 76 };
 77 
 78 class G4DNAEventSet : public IEventSet
 79 {
 80  public:
 81   using Index        = G4VDNAMesh::Index;
 82   using EventSet = std::set<std::unique_ptr<Event>, comparatorEventSet>;
 83   using EventMap = std::unordered_map<Index, EventSet::iterator,G4VDNAMesh::hashFunc>;
 84   G4DNAEventSet();
 85   virtual ~G4DNAEventSet();
 86 
 87   void CreateEvent(const G4double& time, const Index& index,
 88                    Event::ReactionData* pReactionData);
 89   void CreateEvent(const G4double& time, const Index& index,
 90                    std::unique_ptr<Event::JumpingData> jum);
 91 
 92   void AddEvent(std::unique_ptr<Event> pEvent);
 93   void RemoveEventSet()
 94   {
 95     fEventSet.clear();
 96     fEventMap.clear();
 97   }
 98   void RemoveEventOfVoxel(const Index& key);
 99 
100   EventSet::iterator end() { return fEventSet.end(); }
101   EventSet::iterator begin() { return fEventSet.begin(); }
102 
103   EventSet::reverse_iterator rend() { return fEventSet.rend(); }
104   EventSet::reverse_iterator rbegin() { return fEventSet.rbegin(); }
105   EventSet::const_iterator end() const { return fEventSet.end(); }
106   EventSet::const_iterator begin() const { return fEventSet.begin(); }
107   size_t size() { return fEventSet.size(); }
108   G4bool Empty() { return fEventSet.empty(); }
109   void RemoveEvent(EventSet::iterator iter);
110   void PrintEventSet();
111  private:
112   EventSet fEventSet;
113   EventMap fEventMap;
114 };
115 #endif
116