Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEntanglementClipBoard.hh
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//
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
88class G4Track;
89
91
92public:
93
95 : fpParentParticleDefinition(0)
96 , fTrackA(0)
97 , fTrackB(0)
98 , fTrack1Measurement(true)
99 , fTrack2Measurement(true)
100 {}
102
104 {fpParentParticleDefinition = p;}
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
119private:
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
bool G4bool
Definition: G4Types.hh:86
const G4Track * GetTrackB() const
G4ParticleDefinition * GetParentParticleDefinition() const
void SetParentParticleDefinition(G4ParticleDefinition *p)
void SetTrackB(const G4Track *track)
void SetTrackA(const G4Track *track)
const G4Track * GetTrackA() const