14 #include <TClonesArray.h> 15 #include <TDecayChannel.h> 34 #include "Math/GSLMinimizer.h" 35 #include <Math/WrappedParamFunction.h> 37 using namespace genie;
42 Decayer(
"genie::BaryonResonanceDecayer")
48 Decayer(
"genie::BaryonResonanceDecayer", config)
56 for (
unsigned int i = 0; i <
fRParams.size() ; ++i ) {
65 <<
"Running resonance decayer " 69 TObjArrayIter piter(event);
78 int pdg_code = p->
Pdg();
84 if(!this->
ToBeDecayed(pdg_code, status_code))
continue;
87 <<
"Decaying unstable particle: " << p->
Name();
89 bool decayed = this->
Decay(ipos, event);
91 LOG(
"ResonanceDecay",
pWARN) <<
"Resonance not decayed!" ;
92 LOG(
"ResonanceDecay",
pWARN) <<
"Quitting the current event generation thread" ;
97 exception.
SetReason(
"Resonance not decayed");
107 <<
"Done finding & decaying baryon resonances";
117 GHepParticle * decay_particle =
event->Particle(decay_particle_id);
118 if( ! decay_particle) {
120 <<
"Particle to be decayed not in the event record. Particle ud: " << decay_particle_id ;
127 TDecayChannel * selected_decay_channel =
130 if(!selected_decay_channel) {
132 <<
"No decay channel for particle " << decay_particle_id ;
140 bool decayed = this->
DecayExclusive(decay_particle_id, event, selected_decay_channel);
143 delete selected_decay_channel ;
145 if ( ! decayed )
return false ;
149 event->SetWeight(weight);
159 bool & to_be_deleted )
const 162 GHepParticle * decay_particle =
event->Particle(decay_particle_id);
163 if(!decay_particle)
return 0;
166 TLorentzVector decay_particle_p4 = *(decay_particle->
P4());
167 int decay_particle_pdg_code = decay_particle->
Pdg();
170 TParticlePDG * mother =
174 <<
"\n *** The particle with PDG code = " << decay_particle_pdg_code
175 <<
" was not found in PDGLibrary";
179 <<
"Decaying a " << mother->GetName()
184 double W = decay_particle_p4.M();
185 LOG(
"ResonanceDecay",
pINFO) <<
"Available mass W = " <<
W;
188 TObjArray * original_decay_list = mother->DecayList();
190 unsigned int nch = original_decay_list -> GetEntries();
192 << mother->GetName() <<
" has: " << nch <<
" decay channels";
201 TObjArray * actual_decay_list = nullptr ;
203 if ( has_evolved_brs ) {
204 actual_decay_list =
EvolveDeltaBR( decay_particle_pdg_code, original_decay_list, W ) ;
205 if ( ! actual_decay_list )
return nullptr ;
206 nch = actual_decay_list -> GetEntries() ;
207 to_be_deleted = true ;
210 actual_decay_list = original_decay_list ;
211 to_be_deleted = false ;
214 double BR[nch], tot_BR = 0;
216 for(
unsigned int ich = 0; ich < nch; ich++) {
218 TDecayChannel * ch = (TDecayChannel *) actual_decay_list -> At(ich);
224 <<
"Using channel: " << ich
225 <<
" with final state mass = " << fsmass <<
" GeV";
227 tot_BR += ch->BranchingRatio();
231 <<
"Suppresing channel: " << ich
232 <<
" with final state mass = " << fsmass <<
" GeV";
240 <<
"None of the " << nch <<
" decay channels is available @ W = " <<
W;
245 unsigned int ich = 0, sel_ich;
247 double x = tot_BR * rnd->
RndDec().Rndm();
250 }
while (x > BR[ich++]);
252 TDecayChannel * sel_ch = (TDecayChannel *) actual_decay_list -> At(sel_ich);
255 <<
"Selected " << sel_ch->NDaughters() <<
"-particle decay channel (" 256 << sel_ich <<
") has BR = " << sel_ch->BranchingRatio();
258 if ( has_evolved_brs ) {
260 for (
unsigned int i = 0; i < nch; ++i) {
261 if ( sel_ich != i )
delete actual_decay_list -> At(i);
264 delete actual_decay_list ;
274 GHepParticle * decay_particle =
event->Particle(decay_particle_id);
275 if(!decay_particle)
return false ;
278 TLorentzVector decay_particle_p4 = *(decay_particle->
P4());
279 TLorentzVector decay_particle_x4 = *(decay_particle->
X4());
280 int decay_particle_pdg_code = decay_particle->
Pdg();
284 unsigned int nd = ch->NDaughters();
289 for(
unsigned int iparticle = 0; iparticle < nd; iparticle++) {
291 int daughter_code = ch->DaughterPdgCode(iparticle);
295 pdgc[iparticle] = daughter_code;
296 mass[iparticle] = daughter->Mass();
299 <<
"+ daughter[" << iparticle <<
"]: " 300 << daughter->GetName() <<
" (pdg-code = " 301 << pdgc[iparticle] <<
", mass = " << mass[iparticle] <<
")";
316 if ( ! is_permitted )
return false ;
321 for(
int i=0; i<50; i++) {
323 wmax = TMath::Max(wmax,w);
327 <<
"Max phase space gen. weight for current decay: " << wmax;
335 fWeight *= TMath::Max(w/wmax, 1.);
342 bool accept_decay=
false;
351 double gw = wmax * rnd->
RndDec().Rndm();
355 <<
"Current decay weight = " << w <<
" > wmax = " << wmax;
358 <<
"Current decay weight = " << w <<
" / R = " << gw;
360 accept_decay = (gw<=
w);
363 if( accept_decay && is_delta_N_Pi_decay ) {
381 unsigned int pi_id = 0 ;
383 for(
unsigned int iparticle = 0; iparticle < nd; iparticle++) {
407 bool in_nucleus = (target_nucleus!=0);
410 for(
unsigned int id = 0;
id < nd;
id++) {
412 int daughter_pdg_code = pdgc[id];
416 <<
"Adding daughter particle with PDG code = " << pdgc[id]
417 <<
" and mass = " << mass[id] <<
" GeV";
426 daughter_pdg_code, daughter_status_code, decay_particle_id,-1,-1,-1,
427 *daughter_p4, decay_particle_x4);
435 unsigned int nch = decay_list -> GetEntries();
437 std::vector<double> widths( nch, 0. ) ;
440 TDecayChannel *
temp = nullptr ;
442 for (
unsigned int i = 0 ; i < nch ; ++i ) {
444 temp = (TDecayChannel*) decay_list -> At(i) ;
449 if ( tot <= 0. )
return nullptr ;
451 TObjArray *
new_list =
new TObjArray() ;
453 TDecayChannel * update = nullptr ;
455 for (
unsigned int i = 0 ; i < nch ; ++i ) {
457 if ( widths[i] <= 0. ) continue ;
459 temp = (TDecayChannel*) decay_list -> At(i) ;
461 unsigned int nd = temp -> NDaughters() ;
462 std::vector<Int_t>
ds( nd, 0 ) ;
463 for (
unsigned int d = 0 ;
d < nd; ++
d ) {
464 ds[
d] = temp -> DaughterPdgCode(
d) ;
467 update =
new TDecayChannel(
469 temp -> MatrixElementCode(),
475 new_list -> Add( update ) ;
478 new_list -> SetOwner(kFALSE);
502 bool has_pion = false ;
504 int nucleon_id = -1 ;
505 unsigned int nd = ch -> NDaughters() ;
506 for(
unsigned int i = 0 ; i < nd; ++i ) {
551 double m_2 = TMath::Power(m, 2);
554 double mN_2 = TMath::Power( mN, 2);
556 double W_2 = TMath::Power(W, 2);
558 double scaling = 0. ;
563 double m_aux1= TMath::Power( mN + mPion, 2) ;
564 double m_aux2= TMath::Power( mN - mPion, 2) ;
567 double pPi_W = TMath::Sqrt((W_2-m_aux1)*(W_2-m_aux2))/(2*
W);
568 double pPi_m = TMath::Sqrt((m_2-m_aux1)*(m_2-m_aux2))/(2*
m);
570 scaling = TMath::Power( pPi_W / pPi_m , 3 ) ;
576 double Egamma_W = (W_2-mN_2)/(2*W);
577 double Egamma_m = (m_2-mN_2)/(2*m);
580 double fgamma_W = 1./(TMath::Power(1+Egamma_W*Egamma_W/
fFFScaling, 2));
581 double fgamma_m = 1./(TMath::Power(1+Egamma_m*Egamma_m/
fFFScaling, 2));
583 scaling = TMath::Power( Egamma_W / Egamma_m, 3 ) * TMath::Power( fgamma_W / fgamma_m , 2 ) ;
590 return defChWidth * scaling ;
609 GHepParticle * decay_particle =
event->Particle( dec_part_id );
610 TLorentzVector delta_p4 = *(decay_particle->
P4() );
613 TLorentzVector in_lep_p4( * (event -> Probe()-> GetP4()) ) ;
619 TLorentzVector q = in_lep_p4 - out_lep_p4 ;
621 pion.Boost(-delta_p4.BoostVector() );
622 q.Boost(-delta_p4.BoostVector() );
624 TVector3 pion_dir = pion.Vect().Unit() ;
625 TVector3 z_axis = q.Vect().Unit() ;
629 unsigned int q2_index = 0 ;
634 double Q2 = - q.Mag2() ;
644 double c_t = pion_dir*z_axis;
646 w_function = 1. - (
fR33[q2_index] - 0.5)*(3.*c_t*c_t - 1.) ;
654 x[0] = pion_dir.Angle( z_axis ) ;
656 in_lep_p4.Boost(-delta_p4.BoostVector() ) ;
657 out_lep_p4.Boost( -delta_p4.BoostVector() ) ;
660 TVector3 y_axis = in_lep_p4.Vect().Cross( out_lep_p4.Vect() ).Unit() ;
661 TVector3 x_axis = y_axis.Cross(z_axis);
663 TVector3 pion_perp( z_axis.Cross( pion_dir.Cross( z_axis ).Unit() ) ) ;
665 x[1] = pion_perp.Angle( x_axis ) ;
671 if (
fW_max[q2_index] <= 0. ) {
677 return ( aidrnd <= w_function ) ;
691 unsigned int nd = ch->NDaughters();
693 for(
unsigned int iparticle = 0; iparticle < nd; iparticle++) {
695 int daughter_code = ch->DaughterPdgCode(iparticle);
699 double md = daughter->Mass();
703 if(TMath::Abs(daughter_code) == 1114) {
705 <<
"Disabling decay channel containing resonance 1114";;
715 if(!ch)
return false;
717 unsigned int nd = ch->NDaughters();
718 if(nd != 2)
return false;
723 for(
unsigned int iparticle = 0; iparticle < nd; iparticle++) {
725 int daughter_code = ch->DaughterPdgCode(iparticle);
735 if(npi == 1 && nnucl == 1)
return true;
750 <<
"Can decay particle with PDG code = " << pdg_code
751 <<
"? " << ((is_handled)?
"Yes" :
"No");
778 dec_part_pdgc =
abs( dec_part_pdgc ) ;
788 dec_part_pdgc =
abs( dec_part_pdgc ) ;
800 double c_t = TMath::Cos( x[0] ) ;
801 double s_t = TMath::Sin( x[0] ) ;
803 double c_phi = TMath::Cos( x[1 ] );
805 double theta_dep_only = 1. - ( par[0] - 0.5 )*( 3.*c_t*c_t - 1. ) ;
806 double phi_dependency =
kSqrt3 *( 2.*par[1]*s_t*c_t*c_phi + par[2]*s_t*(2.*c_phi*c_phi-1.) ) ;
808 return theta_dep_only - phi_dependency ;
817 ROOT::Math::GSLMinimizer
min( ROOT::Math::kVectorBFGS );
819 min.SetMaxFunctionCalls(1000);
820 min.SetMaxIterations(1000);
823 ROOT::Math::WrappedParamFunction<ROOT::Math::FreeParamMultiFunctionPtr>
f( ( find_maximum ?
829 double step[2] = { 0.00005 *
kPi, 0.00005 * 2 *
kPi } ;
844 }
while( ! min.Minimize() ) ;
846 const double *xs = min.X();
848 double result = find_maximum ? -1 *
f( xs ) :
f( xs ) ;
850 LOG(
"BaryonResonanceDecayer",
pINFO) << (find_maximum ?
"Maximum " :
"Minimum ")
851 <<
"of angular distribution found in ( " 852 << xs[0] <<
", " << xs[1] <<
" ): " 855 LOG(
"BaryonResonanceDecayer",
pDEBUG) <<
"Minimum found in " << min.NCalls() <<
" calls" ;
873 bool invalid_configuration = false ;
879 if (
fR33.size() > 1 ) {
888 invalid_configuration = true ;
889 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"Delta-Q2 and Delta-R33 have wrong sizes" ;
891 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"Delta-R33 -> " <<
fR33.size() ;
897 for (
unsigned int i = 0 ; i <
fR33.size(); ++i ) {
898 if ( (
fR33[i] < -0.5) || (
fR33[i] > 1.) ) {
899 invalid_configuration = true ;
900 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"Delta-R33[" << i <<
"] out of valid range [-0.5 ; 1 ]" ;
901 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"Delta-R33[" << i <<
"] = " <<
fR33[i] ;
908 for (
unsigned int i = 0 ; i <
fR33.size(); ++i ) {
923 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"Delta-R31 or Delta-R3m1 sizes don't match Delta-R33" ;
924 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"R31: " <<
fR31.size()
925 <<
", R3m1: " <<
fR31.size()
926 <<
" while R33: " <<
fR33.size() ;
927 invalid_configuration = true ;
930 for (
unsigned int i = 0; i <
fRParams.size() ; ++i ) {
935 for (
unsigned int i = 0 ; i <
fR33.size(); ++i ) {
942 for (
unsigned int i = 0 ; i <
fR33.size(); ++i ) {
945 if ( temp_min < 0. ) {
946 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"pion angular distribution minimum is negative for Q2 bin " << i ;
947 invalid_configuration = true ;
952 if ( temp_max <= 0. ) {
953 LOG(
"BaryonResonanceDecayer",
pFATAL) <<
"pion angular distribution maximum is non positive for Q2 bin " << i ;
954 invalid_configuration = true ;
964 if ( invalid_configuration ) {
967 <<
"Invalid configuration: Exiting" ;
Baryon resonance decayer module.
virtual bool ToBeDecayed(int pdgc, GHepStatus_t ist) const
const int kPdgP33m1232_DeltaPP
double Weight(void) const
static const double kSqrt3
std::vector< double * > fRParams
bool Decay(int dec_part_id, GHepRecord *event) const
THE MAIN GENIE PROJECT NAMESPACE
static const double kNucleonMass
std::vector< double > fR33
double Q2(const Interaction *const i)
void UnInhibitDecay(int pdgc, TDecayChannel *ch=0) const
static RandomGen * Instance()
Access instance.
const TLorentzVector * P4(void) const
static bool IsDelta(int dec_part_pdgc)
Kinematics * KinePtr(void) const
static const double kPi0Mass
static bool HasEvolvedBRs(int dec_part_pdgc)
int GetParamVect(const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
Handle to load vectors of parameters.
string P4AsString(const TLorentzVector *p)
static const unsigned int kMaxUnweightDecayIterations
bool DecayExclusive(int dec_part_id, GHepRecord *event, TDecayChannel *ch) const
double Mass(Resonance_t res)
resonance mass (GeV)
void Initialize(void) const
double Width(Resonance_t res)
resonance width (GeV)
std::vector< double > fW_max
double EvolveDeltaDecayWidth(int dec_part_pdgc, TDecayChannel *ch, double W) const
enum genie::EResonance Resonance_t
A singleton holding random number generator classes. All random number generation in GENIE should tak...
GHepStatus_t Status(void) const
const int kPdgP33m1232_DeltaP
const int kPdgP33m1232_DeltaM
TDecayChannel * SelectDecayChannel(int dec_part_id, GHepRecord *event, bool &to_be_deleted) const
string Name(void) const
Name that corresponds to the PDG code.
Summary information for an interaction.
double FindDistributionExtrema(unsigned int i, bool find_maximum=false) const
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
const TLorentzVector & FSLeptonP4(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
bool IsPiNDecayChannel(TDecayChannel *ch) const
void ProcessEventRecord(GHepRecord *event) const
bool fRunBefHadroTransp
is invoked before or after FSI?
TGenPhaseSpace fPhaseSpaceGenerator
TObjArray * EvolveDeltaBR(int dec_part_pdgc, TObjArray *decay_list, double W) const
double FinalStateMass(TDecayChannel *ch) const
std::vector< double > fR31
bool fGenerateWeighted
generate weighted or unweighted decays?
Resonance_t FromPdgCode(int pdgc)
PDG code -> resonance id.
Misc GENIE control constants.
std::vector< double > fQ2Thresholds
const int kPdgP33m1232_Delta0
static double MinusPionAngularDist(const double *x, const double *par)
void InhibitDecay(int pdgc, TDecayChannel *ch=0) const
virtual void LoadConfig(void)
void SwitchOnFastForward(void)
static PDGLibrary * Instance(void)
void SetReason(string reason)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
static const double kPionMass
std::vector< double > fR3m1
const TLorentzVector * X4(void) const
void SetStatus(GHepStatus_t s)
static double PionAngularDist(const double *x, const double *par)
Base class for decayer classes. Implements common configuration, allowing users to toggle on/off flag...
bool AcceptPionDecay(TLorentzVector lab_pion, int dec_part_id, const GHepRecord *event) const
virtual ~BaryonResonanceDecayer()
bool IsHandled(int pdgc) const
bool IsBaryonResonance(int pdgc)
is input a baryon resonance?
TParticlePDG * Find(int pdgc, bool must_exist=true)
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
virtual void LoadConfig(void)
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
GENIE's GHEP MC event record.
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
TRandom3 & RndDec(void) const
rnd number generator used by decay models
STDHEP-like event record entry that can fit a particle or a nucleus.
static const double kProtonMass
cet::coded_exception< error, detail::translate > exception
Event finding and building.
enum genie::EGHepStatus GHepStatus_t