Loading [MathJax]/extensions/Safe.js

粒子源定义

用户必须继承抽象基类 G4VUserPrimaryGeneratorAction 来实现粒子源的定义(实现纯虚函数 GeneratePrimaries(G4Eent*))

gun1

建议大家使用 G4ParticleGun ,能够实现灵活度最大化。不建议使用 G4GeneralParticleSource 输入卡的方式。

粒子发射枪

初始化 G4ParticleGun ,以下两行代码通常放在构造函数中。

G4int n_particle = 1;
particleGun = new G4ParticleGun(n_particle);
// 一次发射具有相同属性的 n_particle 个粒子
// 通常设置为 1

粒子定义种类

// 方法一:声明常用粒子
G4ParticleTable* particleTable=G4ParticleTable::GetParticleTable();
G4ParticleDefinition* particle=particleTable->FindParticle("alpha"); 
//“e-”,“e+”,“gamma”,“neutron”,“proton”,…

// 方法二:声明离子
G4IonTable *particleTable=G4IonTable::GetIonTable();
G4ParticleDefinition *particle=particleTable->GetIon(2,4,0.0); //alpha, z,a,ex


particleGun->SetParticleDefinition(particle)

gun2

粒子属性

particleGun->SetParticleCharge(0.);    //电荷态,EM 物理过程自动修正
particleGun->SetParticleEnergy(5.42*MeV);  //动能
particleGun->SetParticleTime(0.); //时间
particleGun->SetParticlePosition(G4ThreeVector(0.,0.,0.)); //初始位置
particleGun->SetParticleMomentumDirection(G4ThreeVector(0.,0.,1.));//动量方向
//particleGun->SetParticleMomentum(G4ParticleMomentum);//动量
//particleGun->SetParticlePolarization(G4ThreeVector);  //极化方向

gun3

粒子发射

particleGun->GeneratePrimaryVertex(anEvent);
// 只有调用了以上代码,才能使该行代码之前的粒子属性设置生效

如果一个事件中发射多个粒子

// ... code 1
// 这里代码省略
particleGun->GeneratePrimaryVertex(anEvent);

// ... code 2
// 这里代码省略
particleGun->GeneratePrimaryVertex(anEvent);

参考代码

void wuPrimaryGeneratorActionAll::GeneratePrimaries(G4Event* anEvent)
{
  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
  G4ParticleDefinition* pp = 0;
  pp = particleTable->FindParticle("gamma");

  if(pp)
    particleGun->SetParticleDefinition(pp);
  else
    G4cout<<"##Null pp in wuPrimaryGeneratorAction::SetParticleGun()"<<G4endl;

  particleGun->SetParticleEnergy(2.0*MeV);  

  // 在 50*50 mm 区域内均匀分布  
  G4double x = (G4UniformRand()-0.5)*50*mm;
  G4double y = (G4UniformRand()-0.5)*50*mm;
  G4double z = 0.;    
  particleGun->SetParticlePosition(G4ThreeVector(x, y, z));  

  // 4pi 各向同性  
  G4double theta= acos(G4UniformRand()-0.5*2);
  G4double phi = G4UniformRand()*2.0*pi;
  G4double cosPX = sin(theta)*cos(phi);
  G4double cosPY = sin(theta)*sin(phi);
  G4double cosPZ = cos(theta);
  G4ThreeVector directPri(cosPX, cosPY, cosPZ);    
  particleGun->SetParticleMomentumDirection(directPri);  

  particleGun->GeneratePrimaryVertex(anEvent);  
}
void wuPrimaryGeneratorActionAll::GeneratePrimaries(G4Event* anEvent)
{
  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
  G4ParticleDefinition* pp = 0;


  pp = particleTable->FindParticle("proton");

  if(pp)
    particleGun->SetParticleDefinition(pp);
  else
    G4cout<<"##Null pp in wuPrimaryGeneratorAction::SetParticleGun()"<<G4endl;

  particleGun->SetParticleEnergy((2.0+8*G4UniformRand())*MeV);  

  // 在 50*50 mm 区域内均匀分布  
  G4double x = (G4UniformRand()-0.5)*50*mm;
  G4double y = (G4UniformRand()-0.5)*50*mm;
  G4double z = 0.;    
  particleGun->SetParticlePosition(G4ThreeVector(x, y, z));  

  // 4pi 各向同性  
  G4double theta= acos(G4UniformRand()-0.5*2);
  G4double phi = G4UniformRand()*2.0*pi;
  G4double cosPX = sin(theta)*cos(phi);
  G4double cosPY = sin(theta)*sin(phi);
  G4double cosPZ = cos(theta);
  G4ThreeVector directPri(cosPX, cosPY, cosPZ);    
  particleGun->SetParticleMomentumDirection(directPri);  

  // 发射第一个粒子  
  particleGun->GeneratePrimaryVertex(anEvent);  

  G4int IonZ = 6;
  G4int IonA = 12;
  G4double IonEstar = 0.0; //exitition energy
  G4double IonQ = 6;
  G4cout<<"ion:Z-A-Q-E*"<<IonZ<<" "<<IonA<<" "<<IonQ<<" "<<IonEstar<<G4endl;
  pp = particleTable->GetIonTable()->GetIon(IonZ, IonA, IonEstar);//4.10.01版本强制 G4IonTable.hh
  // particleGun->SetParticleCharge(IonQ);

  articleGun->SetParticleEnergy((2.0+8*G4UniformRand())*MeV);  
  particleGun->SetParticlePosition(G4ThreeVector(x, y, z));  
  particleGun->SetParticleMomentumDirection(-directPri);    

  // 发射第二个粒子   
  particleGun->GeneratePrimaryVertex(anEvent);  
}