Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/decay/src/G4DecayWithSpin.cc

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 //      GEANT 4 class header file
 28 //
 29 //      History:
 30 //      17 August 2004  P. Gumplinger, T. MacPhail
 31 //      11 April  2008  Kamil Sedlak (PSI), Toni Shiroka (PSI)
 32 // ------------------------------------------------------------
 33 //
 34 #include "G4DecayWithSpin.hh"
 35 
 36 #include "G4Step.hh"
 37 #include "G4Track.hh"
 38 #include "G4DecayTable.hh"
 39 #include "G4MuonDecayChannelWithSpin.hh"
 40 
 41 #include "G4PhysicalConstants.hh"
 42 #include "G4SystemOfUnits.hh"
 43 #include "G4Vector3D.hh"
 44 
 45 #include "G4TransportationManager.hh"
 46 #include "G4PropagatorInField.hh"
 47 #include "G4FieldManager.hh"
 48 #include "G4Field.hh"
 49 
 50 #include "G4Transform3D.hh"
 51 
 52 G4DecayWithSpin::G4DecayWithSpin(const G4String& processName):G4Decay(processName)
 53 {
 54   // set Process Sub Type   
 55   SetProcessSubType(static_cast<int>(DECAY_WithSpin));
 56 
 57 }
 58 
 59 G4DecayWithSpin::~G4DecayWithSpin(){}
 60 
 61 G4VParticleChange* G4DecayWithSpin::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
 62 {
 63   if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
 64        (aTrack.GetTrackStatus() == fStopAndKill )   ){
 65     fParticleChangeForDecay.Initialize(aTrack);
 66     return &fParticleChangeForDecay;
 67   }
 68 
 69 // get particle
 70   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
 71   const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
 72 
 73 // get parent_polarization
 74   G4ThreeVector parent_polarization = aParticle->GetPolarization();
 75 
 76   if(parent_polarization == G4ThreeVector(0,0,0)){
 77     // Generate random polarization direction
 78     G4double cost = 1. - 2.*G4UniformRand();
 79     G4double sint = std::sqrt((1.-cost)*(1.+cost));
 80 
 81     G4double phi = twopi*G4UniformRand();
 82     G4double sinp = std::sin(phi);
 83     G4double cosp = std::cos(phi);
 84 
 85     G4double px = sint*cosp;
 86     G4double py = sint*sinp;
 87     G4double pz = cost;
 88 
 89     parent_polarization.setX(px);
 90     parent_polarization.setY(py);
 91     parent_polarization.setZ(pz);
 92   }
 93 
 94 // decay table
 95   G4DecayTable *decaytable = aParticleDef->GetDecayTable();
 96   if (decaytable != nullptr) {
 97     for (G4int ip=0; ip<decaytable->entries(); ip++){
 98       decaytable->GetDecayChannel(ip)->SetPolarization(parent_polarization);
 99     }
100   }
101 
102   G4ParticleChangeForDecay* pParticleChangeForDecay;
103   pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
104   pParticleChangeForDecay->ProposePolarization(parent_polarization);
105 
106   return pParticleChangeForDecay;
107 }
108   
109 G4VParticleChange* G4DecayWithSpin::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
110 {
111 
112 // get particle
113   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
114   const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
115 
116 // get parent_polarization
117   G4ThreeVector parent_polarization = aParticle->GetPolarization();
118 
119   if(parent_polarization == G4ThreeVector(0,0,0)) {
120     // Generate random polarization direction
121     G4double cost = 1. - 2.*G4UniformRand();
122     G4double sint = std::sqrt((1.-cost)*(1.+cost));
123 
124     G4double phi = twopi*G4UniformRand();
125     G4double sinp = std::sin(phi);
126     G4double cosp = std::cos(phi);
127 
128     G4double px = sint*cosp;
129     G4double py = sint*sinp;
130     G4double pz = cost;
131 
132     parent_polarization.setX(px);
133     parent_polarization.setY(py);
134     parent_polarization.setZ(pz);
135 
136   }else{
137 
138     G4FieldManager* fieldMgr = aStep.GetTrack()->GetVolume()->
139                                      GetLogicalVolume()->GetFieldManager();
140     if (fieldMgr == nullptr) {
141        G4TransportationManager *transportMgr =
142                          G4TransportationManager::GetTransportationManager();
143        G4PropagatorInField* fFieldPropagator = 
144                                         transportMgr->GetPropagatorInField();
145        if (fFieldPropagator) fieldMgr = 
146                                   fFieldPropagator->GetCurrentFieldManager();
147     }
148 
149     const G4Field* field = nullptr;
150     if (fieldMgr != nullptr) field = fieldMgr->GetDetectorField();
151 
152     if ( field != nullptr ) {
153        G4double point[4];
154        point[0] = (aStep.GetPostStepPoint()->GetPosition())[0];
155        point[1] = (aStep.GetPostStepPoint()->GetPosition())[1];
156        point[2] = (aStep.GetPostStepPoint()->GetPosition())[2];
157        point[3] = aTrack.GetGlobalTime();
158 
159        G4double fieldValue[6] ={ 0., 0., 0., 0., 0., 0.};
160        field -> GetFieldValue(point,fieldValue);
161        G4ThreeVector B(fieldValue[0],fieldValue[1],fieldValue[2]);
162 
163        // Call the spin precession only for non-zero mag. field
164        if (B.mag2() > 0.) parent_polarization = 
165                                  Spin_Precession(aStep,B,fRemainderLifeTime);
166 
167     }
168   }
169 
170 // decay table
171   G4DecayTable *decaytable = aParticleDef->GetDecayTable();
172   if ( decaytable != nullptr) {
173     for (G4int ip=0; ip<decaytable->entries(); ip++){
174       decaytable->GetDecayChannel(ip)->SetPolarization(parent_polarization);
175     }
176   } 
177 
178   G4ParticleChangeForDecay* pParticleChangeForDecay;
179   pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
180   pParticleChangeForDecay->ProposePolarization(parent_polarization);
181 
182   return pParticleChangeForDecay;
183 
184 }
185 
186 G4ThreeVector G4DecayWithSpin::Spin_Precession( const G4Step& aStep,
187                                         G4ThreeVector B, G4double deltatime )
188 {
189   G4double Bnorm = std::sqrt(sqr(B[0])  + sqr(B[1]) +sqr(B[2]) );
190 
191   G4double q = aStep.GetTrack()->GetDefinition()->GetPDGCharge();
192   G4double a = 1.165922e-3;
193   G4double s_omega = 8.5062e+7*rad/(s*kilogauss);
194 
195   G4double omega = -(q*s_omega)*(1.+a) * Bnorm;
196 
197   G4double rotationangle = deltatime * omega;
198 
199   G4Transform3D SpinRotation = G4Rotate3D(rotationangle,B.unit());
200 
201   G4Vector3D Spin = aStep.GetTrack() -> GetPolarization();
202 
203   G4Vector3D newSpin = SpinRotation * Spin;
204 
205 #ifdef G4VERBOSE
206   if (GetVerboseLevel()>2) {
207     G4double normspin = std::sqrt(Spin*Spin);
208     G4double normnewspin = std::sqrt(newSpin*newSpin);
209     //G4double cosalpha = Spin*newSpin/normspin/normnewspin;
210     //G4double alpha = std::acos(cosalpha);
211 
212     G4cout << "AT REST::: PARAMETERS " << G4endl;
213     G4cout << "Initial spin  : " << Spin  << G4endl;
214     G4cout << "Delta time    : " << deltatime  << G4endl;
215     G4cout << "Rotation angle: " << rotationangle/rad  << G4endl;
216     G4cout << "New spin      : " << newSpin  << G4endl;
217     G4cout << "Checked norms : " << normspin <<" " << normnewspin << G4endl;
218   }
219 #endif
220 
221   return newSpin;
222 
223 }
224 
225 void G4DecayWithSpin::ProcessDescription(std::ostream& outFile) const
226 {
227   outFile << GetProcessName() 
228     << ": Decay of particles considering parent polarization \n"
229     << "kinematics of daughters are dertermined by DecayChannels \n";
230 }
231