Ignore:
Timestamp:
May 28, 2009, 4:26:57 PM (15 years ago)
Author:
garnier
Message:

maj sur la beta de geant 4.9.3

Location:
trunk/source/processes/hadronic/models/qmd
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/source/processes/hadronic/models/qmd/GNUmakefile

    r819 r1055  
    1 # $Id: GNUmakefile,v 1.1 2007/11/17 16:27:35 dennis Exp $
     1# $Id: GNUmakefile,v 1.2 2008/11/22 06:30:46 tkoi Exp $
    22# -----------------------------------------------------------
    33# GNUmakefile for hadronic library.  Gabriele Cosmo, 18/9/96.
     
    2727            -I$(G4BASE)/processes/hadronic/models/management/include \
    2828            -I$(G4BASE)/processes/hadronic/models/util/include \
     29            -I$(G4BASE)/processes/hadronic/models/theo_high_energy/include \
     30            -I$(G4BASE)/processes/hadronic/models/binary_cascade/include \
    2931            -I$(G4BASE)/processes/hadronic/models/im_r_matrix/include \
     32            -I$(G4BASE)/processes/hadronic/models/chiral_inv_phase_space/interface/include \
     33            -I$(G4BASE)/processes/hadronic/models/chiral_inv_phase_space/body/include \
     34            -I$(G4BASE)/processes/hadronic/models/parton_string/qgsm/include \
     35            -I$(G4BASE)/processes/hadronic/models/parton_string/diffraction/include \
     36            -I$(G4BASE)/processes/hadronic/models/parton_string/hadronization/include \
     37            -I$(G4BASE)/processes/hadronic/models/parton_string/management/include \
    3038            -I$(G4BASE)/processes/hadronic/models/de_excitation/util/include \
    3139            -I$(G4BASE)/processes/hadronic/models/de_excitation/evaporation/include \
  • trunk/source/processes/hadronic/models/qmd/History

    r962 r1055  
    1313     * Please list in reverse chronological order (last date on top)
    1414     ---------------------------------------------------------------
     15
     1631-March-2009 Tatsumi Koi (hadr-qmd-V09-02-01)
     17- Fix bug in gamma(mass zero) partiticpants
     18        src/G4QMDCollision.cc 
     19- Chnage meber object to pointer
     20        include/G4QMDReaction.hh 
     21        src/G4QMDReaction.cc 
     22
     2328-February-2009 Tatsumi Koi (hadr-qmd-V09-02-00)
     24- Fix bug in Absorption
     25        src/G4QMDCollision.cc 
     26- Add Be8 -> Alpha Alpha including opening angle by Q value 
     27        src/G4QMDReaction.cc 
    1528
    162920-November-2008 Tatsumi Koi (hadr-qmd-V09-01-08)
  • trunk/source/processes/hadronic/models/qmd/include/G4QMDReaction.hh

    r962 r1055  
    2828//
    2929//
    30 //      File name:    G4QMDReaction.hh
     30//      File name: G4QMDReaction.hh
    3131//
    3232//      Author: Koi, Tatsumi (tkoi@slac.stanford.edu)       
     
    3737// 081107 Add UnUseGEM (then use the default channel of G4Evaporation)
    3838//            UseFrag (chage criterion of a inelastic reaction)
     39// 090331 Change member shenXS and genspaXS object to pointer
    3940//
    4041
     
    7071      void UseFRAG(){ frag = true; };
    7172
     73      void EnableFermiG(){ heg = true; setHighEnergyModel(); };
     74
    7275   private:
    7376
    7477      void setEvaporationCh();
     78      void setHighEnergyModel();
    7579
    7680      G4QMDMeanField* meanField;
     
    107111      G4double coulomb_collision_pz_targ;
    108112
    109       G4IonsShenCrossSection shenXS;
    110       G4IonsShenCrossSection genspaXS;
     113//090331
     114      G4IonsShenCrossSection* shenXS;
     115      G4GeneralSpaceNNCrossSection* genspaXS;
    111116
    112117      G4bool gem;
    113118      G4bool frag;
     119      G4bool heg;
    114120
    115121};
  • trunk/source/processes/hadronic/models/qmd/src/G4QMDCollision.cc

    r962 r1055  
    2828//        Add several required updating of Mean Filed
    2929//        Modified handling of absorption case by T. Koi
     30// 090126 Fix in absorption case by T. Koi 
     31// 090331 Fix for gamma participant by T. Koi
    3032//
    3133#include "G4QMDCollision.hh"
    32 #include "G4ParticleDefinition.hh"
    3334#include "G4Scatterer.hh"
    3435#include "Randomize.hh"
     
    260261      G4double rmi =  theSystem->GetParticipant( i )->GetMass();
    261262      G4ParticleDefinition* pdi =  theSystem->GetParticipant( i )->GetDefinition();
     263//090331 gamma
     264      if ( pdi->GetPDGMass() == 0.0 ) continue;
    262265
    263266      //std::cout << " p4i00 " << p4i << std::endl;
     
    295298         G4double rmj =  theSystem->GetParticipant( j )->GetMass();
    296299         G4ParticleDefinition* pdj =  theSystem->GetParticipant( j )->GetDefinition();
     300//090331 gamma
     301         if ( pdj->GetPDGMass() == 0.0 ) continue;
    297302
    298303         G4double rr2 = theMeanField->GetRR2( i , j );
     
    332337         G4double pidr = p4i.vect()*dr;
    333338         G4double pjdr = p4j.vect()*dr;
    334  
     339
    335340         G4double aij = 1.0 - ( rmi*rmj /pij ) * ( rmi*rmj /pij );
    336341         G4double bij = pidr / rmi - pjdr*rmi/pij;
     
    660665      else
    661666      {
    662          if ( theMeanField->IsPauliBlocked ( i ) == false )
     667         //if ( theMeanField->IsPauliBlocked ( i ) == false )
     668         //090126                            i-1 cause jth is erased
     669         if ( theMeanField->IsPauliBlocked ( i-1 ) == false )
    663670         {
    664671            //G4cout << "Absorption Happen " << theSystem->GetTotalNumberOfParticipant() << G4endl;
  • trunk/source/processes/hadronic/models/qmd/src/G4QMDReaction.cc

    r962 r1055  
    3131//            UseFrag (chage criterion of a inelastic reaction)
    3232//        Fix bug in nucleon projectiles  by T. Koi   
     33// 090122 Be8 -> Alpha + Alpha
     34// 090331 Change member shenXS and genspaXS object to pointer
    3335//
    3436#include "G4QMDReaction.hh"
     
    4547, frag ( false )
    4648{
     49
     50   //090331
     51   shenXS = new G4IonsShenCrossSection();
     52   //genspaXS = new G4GeneralSpaceNNCrossSection();
    4753   meanField = new G4QMDMeanField();
    4854   collision = new G4QMDCollision();
     
    99105     G4double aTemp = projectile.GetMaterial()->GetTemperature();
    100106
    101      //G4double xs_0 = shenXS.GetCrossSection ( proj_dp , targ_ele , aTemp );
    102      G4double xs_0 = genspaXS.GetCrossSection ( proj_dp , targ_ele , aTemp );
     107     //090331
     108     G4double xs_0 = shenXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
     109     //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
    103110     G4double bmax_0 = std::sqrt( xs_0 / pi );
    104111     //std::cout << "bmax_0 in fm (fermi) " <<  bmax_0/fermi << std::endl;
     
    470477          // Secondary from this nucleus (*it)
    471478          G4ParticleDefinition* pd = (*itt)->GetDefinition();
     479
    472480          G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV );  //in nucleus(*it) rest system
    473481          G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() );  // Back to CM
    474482          G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB 
    475483
    476           G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV ); 
    477           theParticleChange.AddSecondary( dp );
     484
     485//090122
     486          //theParticleChange.AddSecondary( dp );
     487          if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
     488          {
     489             G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV ); 
     490             theParticleChange.AddSecondary( dp );
     491          }
     492          else
     493          {
     494             //Be8 -> Alpha + Alpha + Q
     495             G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
     496             randomized_direction = randomized_direction.unit();
     497             G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
     498             G4double p_decay = std::sqrt ( std::pow(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - std::pow(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
     499             G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 );  //in Be8 rest system
     500             
     501             G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
     502             G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
     503             G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
     504
     505             G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 );  //in Be8 rest system
     506             
     507             G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
     508             G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
     509             G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
     510             
     511             G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV ); 
     512             G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV ); 
     513             theParticleChange.AddSecondary( dp1 );
     514             theParticleChange.AddSecondary( dp2 );
     515          }
     516//090122
    478517
    479518/*
     
    698737
    699738}
    700 
Note: See TracChangeset for help on using the changeset viewer.