Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/management/include/G4VEntanglementClipBoard.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 //
 28 // --------------------------------------------------------------------
 29 // GEANT4 class header file
 30 //
 31 // Class Description:
 32 // Base class for a clipboard for communicating between quantum entangled tracks.
 33 //
 34 // ------------------ G4VEntanglementClipBoard ------------------
 35 //
 36 // Author: J.Allison, May 2017
 37 //
 38 // --------------------------------------------------------------------
 39 
 40 // Usage:
 41 //
 42 // In the method that generates entangled secondaries
 43 // (See, for example, G4eplusAnnihilation::AtRestDoIt)
 44 // Make a shared pointer to a clip board and attach it to the tracks through
 45 // G4EntanglementAuxInfo. That way the clip board lasts the life of both tracks.
 46 // (G4XXXEntanglementClipBoard is your class inherited from this class)
 47 // (See, for example, G4eplusAnnihilationEntanglementClipBoard)
 48 //  auto clipBoard = std::make_shared<G4XXXEntanglementClipBoard>();
 49 // For each secondary
 50 //  G4Track* track = new G4Track(...
 51 //  ...
 52 //  clipBoard->SetTrackA(track);
 53 //  track->SetAuxiliaryTrackInformation(0,new G4EntanglementAuxInfo(clipBoard));
 54 // Then repeat for track B
 55 //
 56 // In the method that does the quantum "measurement"
 57 // (See, for example, G4LivermorePolarizedComptonModel::SampleSecondaries)
 58 //  const auto* auxInfo = fParticleChange->GetCurrentTrack()->GetAuxiliaryTrackInformation(0);
 59 //  if (auxInfo) {
 60 //    const auto* entanglementAuxInfo = dynamic_cast<const G4EntanglementAuxInfo*>(auxInfo);
 61 //    if (entanglementAuxInfo) {
 62 //      auto* clipBoard = dynamic_cast<G4XXXEntanglementClipBoard*>
 63 //      (entanglementAuxInfo->GetEntanglementClipBoard());
 64 //      if (clipBoard) {
 65 //        if (clipBoard->IsTrack1Measurement()) {
 66 //          if (clipBoard->GetTrackB() == fParticleChange->GetCurrentTrack()) {
 67 //            clipBoard->ResetTrack1Measurement();
 68 //            // Store information
 69 //            ...
 70 //          }
 71 //        } else if (clipBoard->IsTrack2Measurement()) {
 72 //          if (clipBoard->GetTrackA() == fParticleChange->GetCurrentTrack()) {
 73 //            clipBoard->ResetTrack2Measurement();
 74 //            // Retrieve information and apply quantum mechanics
 75 //            ...
 76 //          }
 77 //        }
 78 //      }
 79 //    }
 80 //  }
 81 
 82 #ifndef G4VEntanglementClipBoard_hh
 83 #define G4VEntanglementClipBoard_hh
 84 
 85 #include "globals.hh"
 86 
 87 class G4ParticleDefinition;
 88 class G4Track;
 89 
 90 class G4VEntanglementClipBoard {
 91 
 92 public:
 93   
 94   G4VEntanglementClipBoard()
 95   : fpParentParticleDefinition(0)
 96   , fTrackA(0)
 97   , fTrackB(0)
 98   , fTrack1Measurement(true)
 99   , fTrack2Measurement(true)
100   {}
101   virtual ~G4VEntanglementClipBoard() {}
102 
103   void SetParentParticleDefinition(G4ParticleDefinition* p)
104   {fpParentParticleDefinition = p;}
105   G4ParticleDefinition* GetParentParticleDefinition() const
106   {return fpParentParticleDefinition;}
107 
108   void SetTrackA(const G4Track* track) {fTrackA = track;}
109   void SetTrackB(const G4Track* track) {fTrackB = track;}
110   const G4Track* GetTrackA() const {return fTrackA;}
111   const G4Track* GetTrackB() const {return fTrackB;}
112 
113   // The entanglement-sensitive process is responsible for setting this.
114   void ResetTrack1Measurement() {fTrack1Measurement = false;}
115   void ResetTrack2Measurement() {fTrack2Measurement = false;}
116   G4bool IsTrack1Measurement() const {return fTrack1Measurement;}
117   G4bool IsTrack2Measurement() const {return fTrack2Measurement;}
118   
119 private:
120 
121   G4ParticleDefinition* fpParentParticleDefinition;
122 
123   // Use pointer values to identify tracks. We don't know in which order the
124   // tracks will be processed so let us call them A & B.
125   const G4Track* fTrackA;
126   const G4Track* fTrackB;
127 
128   // False until first measurement...
129   G4bool fTrack1Measurement;  // ...of the first encountered track
130   G4bool fTrack2Measurement;  // ...of the second encountered track
131   
132 };
133 
134 #endif
135