musrAtRestSpinRotation Class Reference

#include <musrAtRestSpinRotation.hh>

List of all members.

Public Member Functions

 musrAtRestSpinRotation (const G4String &processName="SpinRotationAtRest")
 musrAtRestSpinRotation (G4VRestProcess &)
virtual ~musrAtRestSpinRotation ()
G4VParticleChange * AtRestDoIt (const G4Track &theTrack, const G4Step &theStep)
void RotateSpin (const G4Step &, G4ThreeVector, G4double)
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *condition)

Static Public Member Functions

static musrAtRestSpinRotationGetInstance ()

Public Attributes

G4double point [4]
G4double B [3]
G4double rotation_angle
G4double itime
G4double ftime
G4double deltatime
G4ThreeVector SpinDirection
G4ParticleChange thePParticleChange
G4VParticleChange theVParticleChange
G4ParticleChangeForTransport theTParticleChange
musrParticleChangeForSR theParticleChange
const G4Field * mfield
G4ThreeVector polar


Detailed Description

Definition at line 14 of file musrAtRestSpinRotation.hh.


Constructor & Destructor Documentation

musrAtRestSpinRotation::musrAtRestSpinRotation const G4String &  processName = "SpinRotationAtRest"  ) 
 

Definition at line 22 of file musrAtRestSpinRotation.cc.

00023   : G4VRestProcess (processName)
00024 {
00025  G4cout << GetProcessName() << " is created "<< G4endl;
00026  pointer = this;
00027 }

musrAtRestSpinRotation::musrAtRestSpinRotation G4VRestProcess &   ) 
 

musrAtRestSpinRotation::~musrAtRestSpinRotation  )  [virtual]
 

Definition at line 30 of file musrAtRestSpinRotation.cc.

00031 {;} 


Member Function Documentation

G4VParticleChange * musrAtRestSpinRotation::AtRestDoIt const G4Track &  theTrack,
const G4Step &  theStep
 

Definition at line 43 of file musrAtRestSpinRotation.cc.

References B, deltatime, ftime, itime, mfield, point, polar, RotateSpin(), and theParticleChange.

00044 {
00045   
00046   theParticleChange.Initialize(theTrack);  
00047   theParticleChange.ProposeTrackStatus(fAlive);
00048   // G4cout<<"AT REST::: globaltime  "<< theTrack.GetGlobalTime()/ns <<G4endl;// getchar();
00049   theParticleChange.ProposePosition(theTrack.GetPosition());
00050   theParticleChange.ProposeMomentumDirection(theTrack.GetMomentumDirection());
00051   theParticleChange.ProposeEnergy(theTrack.GetKineticEnergy());
00052   
00053   // tao :: Get Time
00054   itime = theTrack.GetProperTime();
00055   G4double gtime = theTrack.GetGlobalTime();
00056   ftime = theTrack.GetDynamicParticle()->GetPreAssignedDecayProperTime(); 
00057 
00058 
00059   deltatime = ftime - itime;
00060   theParticleChange.ProposeGlobalTime(deltatime + itime -gtime);
00061   theParticleChange.ProposeProperTime(deltatime);
00062   
00063   polar = aStep.GetTrack()->GetPolarization();
00064 
00065 
00066   
00067   
00068   //  if(aStep.GetTrack()->GetVolume()->GetLogicalVolume()->GetFieldManager())
00069   if(G4TransportationManager::GetTransportationManager()->GetFieldManager())
00070     {
00071       //      G4FieldManager *fMgr = theTrack.GetVolume()->GetLogicalVolume()->GetFieldManager();
00072       G4FieldManager *fMgr=G4TransportationManager::GetTransportationManager()->GetFieldManager();
00073       
00074 
00075 
00076       if(!fMgr->DoesFieldChangeEnergy())//then we have a magnetic field
00077         {
00078           
00079 
00080           //G4cout<<"AT REST::: MAGNETIC FIELD HERE" <<G4endl;// getchar();
00081           
00082           
00083           // tao :: Get Field
00084           point[0]=theTrack.GetPosition().x();
00085           point[1]=theTrack.GetPosition().y();
00086           point[2]=theTrack.GetPosition().z();
00087               
00088           mfield = fMgr->GetDetectorField();       
00089 
00090           mfield->GetFieldValue(point,B);
00091           
00092 
00093           //G4cout <<"AT REST::: MAGNETIC FIELD B="<< B[0] <<" " << B[1] <<" " << B[2] <<G4endl;
00094           
00095 
00096               
00097 #ifdef G4SRVERBOSE
00098           G4cout <<"AT REST::: TIME= proper: "<< itime <<" ;  global: " << gtime  <<" decay: " << ftime <<G4endl;
00099 #endif
00100           
00101           
00102           G4ThreeVector magField(B[0],B[1],B[2]);
00103           RotateSpin(aStep,magField,deltatime);
00104           
00105           
00106 #ifdef G4SRVERBOSE
00107           G4cout<<"AT REST::: spin rotated";
00108 #endif
00109           
00110         } 
00111 
00112 
00113           theParticleChange.ProposePolarization(polar);
00114       
00115 
00116       
00117 #ifdef G4SRVERBOSE
00118       theParticleChange.DumpInfo();
00119 #endif
00120       
00121       
00122       return &theParticleChange;
00123       
00124       // then the AtRestUpdateStep method is called cf LEMuSRParticleChangeForSpinRotation (SR)
00125       
00126     }
00127   
00128   else 
00129     {
00130       //  G4cout <<"No MAGNETIC FIELD do not call AtRestPrecession :: RotateSpin!!! "; 
00131       return &theParticleChange;
00132 
00133 
00134 
00135     }
00136   
00137   
00138 
00139   
00140 }

musrAtRestSpinRotation * musrAtRestSpinRotation::GetInstance  )  [static]
 

Definition at line 35 of file musrAtRestSpinRotation.cc.

00036 {
00037   return pointer;
00038 }

G4double musrAtRestSpinRotation::GetMeanLifeTime const G4Track &  ,
G4ForceCondition *  condition
[inline]
 

Definition at line 45 of file musrAtRestSpinRotation.hh.

00046      {
00047        *condition = Forced;
00048        return DBL_MAX;
00049      }

void musrAtRestSpinRotation::RotateSpin const G4Step &  ,
G4ThreeVector  ,
G4double 
 

Definition at line 144 of file musrAtRestSpinRotation.cc.

References polar, and rotation_angle.

Referenced by AtRestDoIt().

00145 {
00146   
00147   
00148   G4Transform3D Spin_rotation;
00149  
00150   G4double Bnorm = sqrt(sqr(B[0])  + sqr(B[1]) +sqr(B[2]) );
00151 
00152 
00153   //  G4cout<< "AT REST::: PARAMETERS\n" 
00154   //    << "Magnetic Field Norm  : " << G4BestUnit(Bnorm,"Magnetic flux density") <<"\n";
00155     
00156   G4double omega,q,a,fqz;
00157   
00158   q= aStep.GetTrack()->GetDefinition()->GetPDGCharge();
00159   a= 1.165922e-3;
00160   fqz = 8.5062e+7*rad/(s*kilogauss);
00161   
00162 #ifdef G4SRVERBOSE
00163   G4cout<< "AT REST::: PARAMETERS\n" 
00164         << "Charge  : " << q <<"\n";
00165 #endif
00166   
00167   omega= - (fqz)*(1.+a) * Bnorm;
00168 
00169 #ifdef G4SRVERBOSE
00170   G4cout<< "AT REST::: PARAMETERS\n" 
00171         << "Frequency: " << G4BestUnit(omega/(2*M_PI*rad),"Frequency") <<"\n";
00172 #endif
00173    
00174   
00175   rotation_angle = deltatime*omega; 
00176   
00177 
00178   Spin_rotation= G4Rotate3D(rotation_angle,B/Bnorm);
00179   
00180   HepVector3D spin = aStep.GetTrack()->GetPolarization();
00181   HepVector3D newspin;
00182   newspin = Spin_rotation*spin;
00183   
00184   G4double x,y,z,alpha;
00185   x = sqrt(spin*spin);
00186   y = sqrt(newspin*newspin);
00187   z = spin*newspin/x/y;
00188   alpha = acos(z);
00189   
00190   #ifdef G4SRVERBOSE
00191     G4cout<< "AT REST::: PARAMETERS\n" 
00192   << "Initial spin  : " << spin <<"\n"
00193   << "Delta time    : " << deltatime <<"\n"
00194   << "Rotation angle: " << rotation_angle/rad <<"\n"
00195   << "New spin      : " << newspin <<"\n"
00196   << "Checked norms : " << x <<" " <<y <<" \n"
00197   << G4endl; 
00198   #endif
00199   
00200   polar = newspin;
00201   
00202 }


Member Data Documentation

G4double musrAtRestSpinRotation::B[3]
 

Definition at line 33 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4double musrAtRestSpinRotation::deltatime
 

Definition at line 35 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4double musrAtRestSpinRotation::ftime
 

Definition at line 35 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4double musrAtRestSpinRotation::itime
 

Definition at line 35 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

const G4Field* musrAtRestSpinRotation::mfield
 

Definition at line 43 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4double musrAtRestSpinRotation::point[4]
 

Definition at line 32 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4ThreeVector musrAtRestSpinRotation::polar
 

Definition at line 50 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt(), and RotateSpin().

G4double musrAtRestSpinRotation::rotation_angle
 

Definition at line 33 of file musrAtRestSpinRotation.hh.

Referenced by RotateSpin().

G4ThreeVector musrAtRestSpinRotation::SpinDirection
 

Definition at line 36 of file musrAtRestSpinRotation.hh.

musrParticleChangeForSR musrAtRestSpinRotation::theParticleChange
 

Definition at line 42 of file musrAtRestSpinRotation.hh.

Referenced by AtRestDoIt().

G4ParticleChange musrAtRestSpinRotation::thePParticleChange
 

Definition at line 39 of file musrAtRestSpinRotation.hh.

G4ParticleChangeForTransport musrAtRestSpinRotation::theTParticleChange
 

Definition at line 41 of file musrAtRestSpinRotation.hh.

G4VParticleChange musrAtRestSpinRotation::theVParticleChange
 

Definition at line 40 of file musrAtRestSpinRotation.hh.


The documentation for this class was generated from the following files:
Generated on Mon Mar 27 12:19:54 2006 for MUSR by  doxygen 1.4.6