G4NeutronHPChannel Class Reference

#include <G4NeutronHPChannel.hh>


Public Member Functions

 G4NeutronHPChannel ()
 ~G4NeutronHPChannel ()
G4double GetXsec (G4double energy)
G4double GetWeightedXsec (G4double energy, G4int isoNumber)
G4double GetFSCrossSection (G4double energy, G4int isoNumber)
G4bool IsActive (G4int isoNumber)
G4bool HasFSData (G4int isoNumber)
G4bool HasAnyData (G4int isoNumber)
G4bool Register (G4NeutronHPFinalState *theFS)
void Init (G4Element *theElement, const G4String dirName)
void Init (G4Element *theElement, const G4String dirName, const G4String fsType)
void UpdateData (G4int A, G4int Z, G4int index, G4double abundance)
void UpdateData (G4int A, G4int Z, G4int M, G4int index, G4double abundance)
void Harmonise (G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack, G4int isoNumber=-1)
G4int GetNiso ()
G4double GetN (G4int i)
G4double GetZ (G4int i)
G4bool HasDataInAnyFinalState ()


Detailed Description

Definition at line 54 of file G4NeutronHPChannel.hh.


Constructor & Destructor Documentation

G4NeutronHPChannel::G4NeutronHPChannel (  )  [inline]

Definition at line 58 of file G4NeutronHPChannel.hh.

00059   {
00060     theChannelData = new G4NeutronHPVector; 
00061     theBuffer = 0;
00062     theIsotopeWiseData = 0;
00063     theFinalStates = 0;
00064     active = 0;
00065     registerCount = -1;
00066   }

G4NeutronHPChannel::~G4NeutronHPChannel (  )  [inline]

Definition at line 68 of file G4NeutronHPChannel.hh.

00069   {
00070     delete theChannelData; 
00071     // Following statement disabled to avoid SEGV
00072     // theBuffer is also deleted as "theChannelData" in
00073     // ~G4NeutronHPIsoData.  FWJ 06-Jul-1999
00074     //if(theBuffer != 0) delete theBuffer; 
00075     if(theIsotopeWiseData != 0) delete [] theIsotopeWiseData;
00076     // Deletion of FinalStates disabled to avoid endless looping
00077     // in the destructor heirarchy.  FWJ 06-Jul-1999
00078     //if(theFinalStates != 0)
00079     //{
00080     //  for(i=0; i<niso; i++)
00081     //  {
00082     //    delete theFinalStates[i];
00083     //  }
00084     //  delete [] theFinalStates;
00085     //}
00086     // FWJ experiment
00087     //if(active!=0) delete [] active;
00088 // T.K. 
00089    if ( theFinalStates != 0 )
00090    {
00091       for ( G4int i = 0 ; i < niso ; i++ )
00092       {
00093          delete theFinalStates[i];
00094       }
00095       delete [] theFinalStates;
00096    }
00097    if ( active != 0 ) delete [] active;
00098     
00099   }


Member Function Documentation

G4HadFinalState * G4NeutronHPChannel::ApplyYourself ( const G4HadProjectile theTrack,
G4int  isoNumber = -1 
)

Definition at line 211 of file G4NeutronHPChannel.cc.

References G4NeutronHPFinalState::ApplyYourself(), G4UniformRand, G4NeutronHPManager::GetInstance(), G4NeutronHPFinalState::GetN(), G4NeutronHPManager::GetReactionWhiteBoard(), G4NeutronHPThermalBoost::GetThermalEnergy(), GetXsec(), G4NeutronHPFinalState::GetZ(), HasAnyData(), G4NeutronHPReactionWhiteBoard::SetTargA(), and G4NeutronHPReactionWhiteBoard::SetTargZ().

Referenced by G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPChannelList::ApplyYourself(), and G4FissLib::ApplyYourself().

00212   {
00213 //    G4cout << "G4NeutronHPChannel::ApplyYourself+"<<niso<<G4endl;
00214     if ( anIsotope != -1 ) 
00215     {
00216        //Inelastic Case
00217        //G4cout << "G4NeutronHPChannel Inelastic Case" 
00218        //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
00219        G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->SetTargA( (G4int)this->GetN(anIsotope) ); 
00220        G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->SetTargZ( (G4int)this->GetZ(anIsotope) ); 
00221        return theFinalStates[anIsotope]->ApplyYourself(theTrack);
00222     }
00223     G4double sum=0;
00224     G4int it=0;
00225     G4double * xsec = new G4double[niso];
00226     G4NeutronHPThermalBoost aThermalE;
00227     for (G4int i=0; i<niso; i++)
00228     {
00229       if(theFinalStates[i]->HasAnyData())
00230       {
00231         xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
00232                                                                            theFinalStates[i]->GetN(),
00233                                                                            theFinalStates[i]->GetZ(),
00234                                                                            theTrack.GetMaterial()->GetTemperature()));
00235         sum += xsec[i];
00236       }
00237       else
00238       {
00239         xsec[i]=0;
00240       }
00241     } 
00242     if(sum == 0) 
00243     {
00244 //      G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
00245 //      G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
00246       it = static_cast<G4int>(niso*G4UniformRand());
00247     }
00248     else
00249     {
00250 //      G4cout << "Are we still here? "<<sum<<G4endl;
00251 //      G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
00252       G4double random = G4UniformRand();
00253       G4double running=0;
00254 //      G4cout << "G4NeutronHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
00255 //      G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
00256       for (G4int ix=0; ix<niso; ix++)
00257       {
00258         running += xsec[ix];
00259         //if(random<=running/sum) 
00260         if( sum == 0 || random <= running/sum ) 
00261         { 
00262           it = ix;
00263           break;
00264         }
00265       }
00266       if(it==niso) it--;
00267     }
00268     delete [] xsec;
00269     G4HadFinalState * theFinalState=0;
00270     while(theFinalState==0)
00271     {
00272 //    G4cout << "TESTHP 24 it="<<it<<G4endl;
00273       theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
00274     }
00275 //    G4cout <<"THE IMPORTANT RETURN"<<G4endl;
00276       //G4cout << "TK G4NeutronHPChannel Elastic, Capture and Fission Cases " 
00277       //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
00278     G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->SetTargA( (G4int)this->GetN(it) ); 
00279     G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->SetTargZ( (G4int)this->GetZ(it) ); 
00280     return theFinalState;
00281   }

G4double G4NeutronHPChannel::GetFSCrossSection ( G4double  energy,
G4int  isoNumber 
)

Definition at line 53 of file G4NeutronHPChannel.cc.

References G4NeutronHPFinalState::GetXsec().

00054   {
00055     return theFinalStates[isoNumber]->GetXsec(energy);
00056   }

G4double G4NeutronHPChannel::GetN ( G4int  i  )  [inline]

Definition at line 129 of file G4NeutronHPChannel.hh.

References G4NeutronHPFinalState::GetN().

Referenced by G4NeutronHPChannelList::ApplyYourself().

00129 {return theFinalStates[i]->GetN();}

G4int G4NeutronHPChannel::GetNiso (  )  [inline]

Definition at line 127 of file G4NeutronHPChannel.hh.

Referenced by G4NeutronHPorLFission::GetNiso(), G4NeutronHPorLElastic::GetNiso(), and G4NeutronHPorLCapture::GetNiso().

00127 {return niso;}

G4double G4NeutronHPChannel::GetWeightedXsec ( G4double  energy,
G4int  isoNumber 
)

Definition at line 48 of file G4NeutronHPChannel.cc.

References G4NeutronHPIsoData::GetXsec().

00049   {
00050     return theIsotopeWiseData[isoNumber].GetXsec(energy);
00051   }

G4double G4NeutronHPChannel::GetXsec ( G4double  energy  ) 

Definition at line 43 of file G4NeutronHPChannel.cc.

References G4NeutronHPVector::GetXsec().

Referenced by G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), ApplyYourself(), G4FissLib::ApplyYourself(), G4NeutronHPorLFissionData::GetCrossSection(), G4NeutronHPorLElasticData::GetCrossSection(), and G4NeutronHPorLCaptureData::GetCrossSection().

00044   {
00045     return std::max(0., theChannelData->GetXsec(energy));
00046   }

G4double G4NeutronHPChannel::GetZ ( G4int  i  )  [inline]

Definition at line 130 of file G4NeutronHPChannel.hh.

References G4NeutronHPFinalState::GetZ().

Referenced by G4NeutronHPChannelList::ApplyYourself(), and Register().

00130 {return theFinalStates[i]->GetZ();}

void G4NeutronHPChannel::Harmonise ( G4NeutronHPVector *&  theStore,
G4NeutronHPVector theNew 
)

Definition at line 166 of file G4NeutronHPChannel.cc.

References G4NeutronHPVector::GetEnergy(), G4NeutronHPVector::GetVectorLength(), G4NeutronHPVector::GetXsec(), CLHEP::detail::n, and G4NeutronHPVector::SetData().

Referenced by UpdateData().

00167   {
00168     G4int s_tmp = 0, n=0, m_tmp=0;
00169     G4NeutronHPVector * theMerge = new G4NeutronHPVector;
00170     G4NeutronHPVector * anActive = theStore;
00171     G4NeutronHPVector * aPassive = theNew;
00172     G4NeutronHPVector * tmp;
00173     G4int a = s_tmp, p = n, t;
00174     while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength())
00175     {
00176       if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
00177       {
00178         G4double xa  = anActive->GetEnergy(a);
00179         theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
00180         m_tmp++;
00181         a++;
00182         G4double xp = aPassive->GetEnergy(p);
00183         if( std::abs(std::abs(xp-xa)/xa)<0.001 )
00184         {
00185           p++;
00186         }
00187       } else {
00188         tmp = anActive; t=a;
00189         anActive = aPassive; a=p;
00190         aPassive = tmp; p=t;
00191       }
00192     }
00193     while (a!=anActive->GetVectorLength())
00194     {
00195       theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
00196       a++;
00197     }
00198     while (p!=aPassive->GetVectorLength())
00199     {
00200       if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
00201         theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
00202       p++;
00203     }
00204     delete theStore;
00205     theStore = theMerge;
00206   }

G4bool G4NeutronHPChannel::HasAnyData ( G4int  isoNumber  )  [inline]

Definition at line 111 of file G4NeutronHPChannel.hh.

References G4NeutronHPFinalState::HasAnyData().

Referenced by ApplyYourself(), HasDataInAnyFinalState(), and UpdateData().

00111 { return theFinalStates[isoNumber]->HasAnyData(); }

G4bool G4NeutronHPChannel::HasDataInAnyFinalState (  )  [inline]

Definition at line 132 of file G4NeutronHPChannel.hh.

References HasAnyData().

Referenced by Register().

00133   {
00134     G4bool result = false;
00135     G4int i;
00136     for(i=0; i<niso; i++)
00137     {
00138       if(theFinalStates[i]->HasAnyData()) result = true;
00139     }
00140     return result;
00141   }

G4bool G4NeutronHPChannel::HasFSData ( G4int  isoNumber  )  [inline]

Definition at line 109 of file G4NeutronHPChannel.hh.

References G4NeutronHPFinalState::HasFSData().

00109 { return theFinalStates[isoNumber]->HasFSData(); }

void G4NeutronHPChannel::Init ( G4Element theElement,
const G4String  dirName,
const G4String  fsType 
)

Definition at line 59 of file G4NeutronHPChannel.cc.

References Init().

00060   {
00061     theFSType = aFSType;
00062     Init(anElement, dirName);
00063   }

void G4NeutronHPChannel::Init ( G4Element theElement,
const G4String  dirName 
)

Definition at line 65 of file G4NeutronHPChannel.cc.

Referenced by Init(), and G4NeutronHPChannelList::Register().

00066   {
00067     theDir = dirName;
00068     theElement = anElement;
00069   }

G4bool G4NeutronHPChannel::IsActive ( G4int  isoNumber  )  [inline]

Definition at line 107 of file G4NeutronHPChannel.hh.

00107 { return active[isoNumber]; }

G4bool G4NeutronHPChannel::Register ( G4NeutronHPFinalState theFS  ) 

Definition at line 71 of file G4NeutronHPChannel.cc.

References G4lrint(), G4StableIsotopes::GetAbundance(), G4StableIsotopes::GetFirstIsotope(), G4Element::GetIsotope(), G4StableIsotopes::GetIsotopeNucleonCount(), G4Isotope::Getm(), G4Isotope::GetN(), G4StableIsotopes::GetNumberOfIsotopes(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), GetZ(), G4Element::GetZ(), HasDataInAnyFinalState(), G4NeutronHPFinalState::New(), G4NeutronHPFinalState::SetA_Z(), and UpdateData().

Referenced by G4FissLib::G4FissLib(), and G4NeutronHPChannelList::Register().

00072   {
00073     registerCount++;
00074     G4int Z = G4lrint(theElement->GetZ());
00075 
00076     Z = Z-registerCount;
00077     if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case 
00078     if ( Z < 1 ) return false; 
00079 /*
00080     if(registerCount<5)
00081     {
00082       Z = Z-registerCount;
00083     }
00084 */
00085     //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
00086     // Bug fix by TK on behalf of AH
00087     if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
00088     G4int count = 0;
00089     if(registerCount==0) count = theElement->GetNumberOfIsotopes();
00090     if(count == 0||registerCount!=0) count +=
00091          theStableOnes.GetNumberOfIsotopes(Z);
00092     niso = count;
00093     delete [] theIsotopeWiseData;
00094     theIsotopeWiseData = new G4NeutronHPIsoData [niso];
00095     delete [] active;
00096     active = new G4bool[niso];
00097 
00098     delete [] theFinalStates;
00099     theFinalStates = new G4NeutronHPFinalState * [niso];
00100     delete theChannelData;
00101     theChannelData = new G4NeutronHPVector; 
00102     for(G4int i=0; i<niso; i++)
00103     {
00104       theFinalStates[i] = theFS->New();
00105     }
00106     count = 0;
00107     G4int nIsos = niso;
00108     if(theElement->GetNumberOfIsotopes()!=0&&registerCount==0)
00109     {
00110       for (G4int i1=0; i1<nIsos; i1++)
00111       {
00112         // G4cout <<" Init: normal case"<<G4endl;
00113         G4int A = theElement->GetIsotope(i1)->GetN();
00114         G4int M = theElement->GetIsotope(i1)->Getm();
00115         G4double frac = theElement->GetRelativeAbundanceVector()[i1]/perCent;
00116         //theFinalStates[i1]->SetA_Z(A, Z);
00117         //UpdateData(A, Z, count++, frac);
00118         theFinalStates[i1]->SetA_Z(A, Z, M);
00119         UpdateData(A, Z, M, count++, frac);
00120       }
00121     } else {
00122       //G4cout <<" Init: mean case: "
00123       //       <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
00124         //     <<Z<<" "<<theElement
00125         //     << G4endl;
00126       G4int first = theStableOnes.GetFirstIsotope(Z);
00127       for(G4int i1=0; 
00128         i1<theStableOnes.GetNumberOfIsotopes(Z);
00129         i1++)
00130       {
00131         G4int A = theStableOnes.GetIsotopeNucleonCount(first+i1);
00132         G4double frac = theStableOnes.GetAbundance(first+i1);
00133         theFinalStates[i1]->SetA_Z(A, Z);
00134         UpdateData(A, Z, count++, frac);
00135       }
00136     }
00137     G4bool result = HasDataInAnyFinalState();
00138     return result;
00139   }

void G4NeutronHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  M,
G4int  index,
G4double  abundance 
)

Definition at line 142 of file G4NeutronHPChannel.cc.

References G4NeutronHPIsoData::FillChannelData(), G4NeutronHPFinalState::GetXsec(), Harmonise(), HasAnyData(), G4NeutronHPIsoData::Init(), G4NeutronHPFinalState::Init(), G4NeutronHPIsoData::MakeChannelData(), and G4NeutronHPVector::Times().

00143   {
00144     //theFinalStates[index]->Init(A, Z, theDir, theFSType);
00145     theFinalStates[index]->Init(A, Z, M, theDir, theFSType);
00146     if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
00147 
00148     // the above has put the X-sec into the FS
00149     theBuffer = 0;
00150     if(theFinalStates[index]->HasXsec())
00151     {
00152       theBuffer = theFinalStates[index]->GetXsec();
00153       theBuffer->Times(abundance/100.);
00154       theIsotopeWiseData[index].FillChannelData(theBuffer);
00155     }
00156     else // get data from CrossSection directory
00157     {
00158       G4String tString = "/CrossSection";
00159       //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
00160       active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
00161       if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
00162     }
00163     if(theBuffer != 0) Harmonise(theChannelData, theBuffer);
00164   }

void G4NeutronHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  index,
G4double  abundance 
) [inline]

Definition at line 120 of file G4NeutronHPChannel.hh.

Referenced by Register().

00120 { G4int M = 0; UpdateData( A, Z, M, index, abundance); };


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:36 2013 for Geant4 by  doxygen 1.4.7