diff --git a/inc/LauAbsResonance.hh b/inc/LauAbsResonance.hh index cc16996..335ee87 100644 --- a/inc/LauAbsResonance.hh +++ b/inc/LauAbsResonance.hh @@ -1,578 +1,579 @@ /* Copyright 2004 University of Warwick Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except in compliance with the License. You may obtain a copy of the License at http://www.apache.org/licenses/LICENSE-2.0 Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License. */ /* Laura++ package authors: John Back Paul Harrison Thomas Latham */ /*! \file LauAbsResonance.hh \brief File containing declaration of LauAbsResonance class. */ /*! \class LauAbsResonance \brief Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc.) Abstract Class for defining the type for all classes used to model resonances in the Dalitz plot, such as Breit-Wigner functions. In addition, some common functionality is implemented, including data such as the mass and width of the desired state. */ #ifndef LAU_ABS_RESONANCE #define LAU_ABS_RESONANCE #include "TString.h" #include "LauBlattWeisskopfFactor.hh" #include "LauComplex.hh" #include "LauParameter.hh" class LauDaughters; class LauKinematics; class LauResonanceInfo; class LauAbsResonance { public: //! Define the allowed resonance types enum LauResonanceModel { BW, /*!< simple Breit-Wigner */ RelBW, /*!< relativistic Breit-Wigner */ GS, /*!< a modified Breit-Wigner from Gounaris-Sakurai */ Flatte, /*!< Flatte or coupled-channel Breit-Wigner */ Sigma, /*!< special shape for the sigma or f_0(600) */ Kappa, /*!< special shape for the kappa, a low-mass Kpi scalar */ Dabba, /*!< special shape for the dabba, a low-mass Dpi scalar */ LASS, /*!< the LASS amplitude to describe the Kpi S-wave */ LASS_BW, /*!< the resonant part of the LASS amplitude */ LASS_NR, /*!< the nonresonant part of the LASS amplitude */ EFKLLM, /*!< a form-factor-based description of the Kpi S-wave */ KMatrix, /*!< S-wave description using K-matrix and P-vector */ FlatNR, /*!< a uniform nonresonant amplitude */ NRModel, /*!< a theoretical model nonresonant amplitude */ BelleNR, /*!< an empirical exponential nonresonant amplitude */ + RescatterThreshold, /*!< a theoretical model to account for rescattering from a nearby-threshold channel */ PowerLawNR, /*!< an empirical power law nonresonant amplitude */ BelleSymNR, /*!< an empirical exponential nonresonant amplitude for symmetrised DPs */ BelleSymNRNoInter, /*!< an empirical exponential nonresonant amplitude for symmetrised DPs without interference */ TaylorNR, /*!< an empirical Taylor expansion nonresonant amplitude for symmetrised DPs */ PolNR, /*!< an empirical polynomial nonresonant amplitude */ Pole, /*!< scalar Pole lineshape */ PolarFFNR, /*!< Polar Form Factor nonresonant amplitude */ PolarFFSymNR, /*!< Polar Form Factor nonresonant amplitude for symmetrised DPs */ PolarFFSymNRNoInter, /*!< Polar Form Factor nonresonant amplitude for symmetrised DPs without interference */ Rescattering, /*!< KK-PiPi inelastic scattering amplitude */ Rescattering2, /*!< KK-PiPi inelastic scattering amplitude */ RescatteringNoInter, /*!< KK-PiPi inelastic scattering amplitude */ MIPW_MagPhase, /*!< a model independent partial wave - magnitude and phase representation */ MIPW_RealImag, /*!< a model independent partial wave - real and imaginary part representation */ GaussIncoh, /*!< an incoherent Gaussian shape */ RhoOmegaMix_GS, /*!< mass mixing model using GS for res 1 and RBW for res 2 */ RhoOmegaMix_RBW, /*!< mass mixing model using two RBWs */ RhoOmegaMix_GS_1, /*!< mass mixing model using GS for res 1 and RBW for res 2, with denominator factor = 1 */ RhoOmegaMix_RBW_1 /*!< mass mixing model using two RBWs, with denominator factor = 1 */ }; //! Define the allowed spin formalisms enum LauSpinType { Zemach_P, /*!< Zemach tensor formalism, bachelor momentum in resonance rest frame */ Zemach_Pstar, /*!< Zemach tensor formalism, bachelor momentum in parent rest frame */ Covariant, /*!< Covariant tensor formalism */ Legendre /*!< Legendre polynomials only */ }; //! Is the resonance model incoherent? /*! \param [in] model the resonance model \return true if the model is incoherent */ static bool isIncoherentModel(LauResonanceModel model); //! Constructor (for use by standard resonances) /*! \param [in] resInfo the object containing information on the resonance name, mass, width, spin, charge, etc. \param [in] resPairAmpInt the number of the daughter not produced by the resonance \param [in] daughters the daughter particles */ LauAbsResonance(LauResonanceInfo* resInfo, const Int_t resPairAmpInt, const LauDaughters* daughters); //! Constructor (for use by K-matrix components) /*! \param [in] resName the name of the component \param [in] resPairAmpInt the number of the daughter not produced by the resonance \param [in] daughters the daughter particles */ LauAbsResonance(const TString& resName, const Int_t resPairAmpInt, const LauDaughters* daughters); //! Destructor virtual ~LauAbsResonance(); //! Initialise the model virtual void initialise() = 0; //! Calculate the complex amplitude /*! \param [in] kinematics the kinematic variables of the current event \return the complex amplitude */ virtual LauComplex amplitude(const LauKinematics* kinematics); //! Get the resonance model type /*! \return the resonance model type */ virtual LauResonanceModel getResonanceModel() const = 0; //! Get the spin type /*! \return the spin formalism */ LauSpinType getSpinType() const {return spinType_;} //! Get the name of the resonance /*! \return the resonance name */ const TString& getResonanceName() const {return resName_;} //! Get the name of the resonance /*! \return the resonance name */ const TString& getSanitisedName() const {return sanitisedName_;} //! Get the integer to identify which DP axis the resonance belongs to /*! \return the DP axis identification number, the ID of the bachelor */ Int_t getPairInt() const {return resPairAmpInt_;} //! Get the spin of the resonance /*! \return the resonance spin */ Int_t getSpin() const {return resSpin_;} //! Get the charge of the resonance /*! \return the resonance charge */ Int_t getCharge() const {return resCharge_;} //! Get the mass of the resonance /*! \return the resonance mass */ Double_t getMass() const {return (resMass_!=0) ? resMass_->unblindValue() : -1.0;} //! Get the width of the resonance /*! \return the resonance width */ Double_t getWidth() const {return (resWidth_!=0) ? resWidth_->unblindValue() : -1.0;} //! Get the mass parameter of the resonance /*! \return the resonance mass parameter */ LauParameter* getMassPar() {return resMass_;} //! Get the width parameter of the resonance /*! \return the resonance width parameter */ LauParameter* getWidthPar() {return resWidth_;} //! Retrieve the resonance parameters, e.g. so that they can be loaded into a fit /*! \return floating parameters of the resonance */ virtual const std::vector& getFloatingParameters() { return this->getParameters(); }; //! Is the amplitude pre-symmetrised? /*! The default value is kFALSE, so pre-symmetrised lineshapes should override this. \return whether the amplitude is already symmetrised */ virtual Bool_t preSymmetrised() const {return kFALSE;} //! Get the helicity flip flag /*! \return the flip helicity flag */ Bool_t flipHelicity() const {return flipHelicity_;} //! Set the helicity flip flag /*! \param [in] boolean the helicity flip status */ void flipHelicity(const Bool_t boolean) {flipHelicity_ = boolean;} //! Get the ignore momenta flag /*! Whether to ignore the momentum factors in both the spin factor and the mass-dependent width \return the ignore momenta flag */ Bool_t ignoreMomenta() const {return ignoreMomenta_;} //! Set the ignore momenta flag /*! Whether to ignore the momentum factors in both the spin factor and the mass-dependent width \param [in] boolean the ignore momenta status */ void ignoreMomenta(const Bool_t boolean) {ignoreMomenta_ = boolean;} //! Get the ignore spin flag /*! Whether to set the spinTerm to unity always \return the ignore spin flag */ Bool_t ignoreSpin() const {return ignoreSpin_;} //! Set the ignore spin flag /*! Whether to set the spinTerm to unity always \param [in] boolean the ignore spin status */ void ignoreSpin(const Bool_t boolean) {ignoreSpin_ = boolean;} //! Get the ignore barrier factor scaling flag /*! Whether to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width \return the ignore barrier amplitude scaling flag */ Bool_t ignoreBarrierScaling() const {return ignoreBarrierScaling_;} //! Set the ignore barrier factor scaling flag /*! Whether to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width \param [in] boolean the ignore barrier factor scaling status */ void ignoreBarrierScaling(const Bool_t boolean) {ignoreBarrierScaling_ = boolean;} //! Allow the mass, width and spin of the resonance to be changed /*! Negative values wil be ignored, so if, for example, you want to only change the spin you can provide negative values for the mass and width \param [in] newMass new value of the resonance mass \param [in] newWidth new value of the resonance width \param [in] newSpin new value of the resonance spin */ void changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin); //! Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed /*! Negative values wil be ignored, so if, for example, you want to only change the parent radius you can provide a negative value for the resonance radius \param [in] resRadius new value of the resonance radius \param [in] parRadius new value of the parent radius */ void changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius); //! Set value of the various parameters /*! \param [in] name the name of the parameter to be changed \param [in] value the new parameter value */ virtual void setResonanceParameter(const TString& name, const Double_t value); //! Allow the various parameters to float in the fit /*! \param [in] name the name of the parameter to be floated */ virtual void floatResonanceParameter(const TString& name); //! Access the given resonance parameter /*! \param [in] name the name of the parameter \return the corresponding parameter */ virtual LauParameter* getResonanceParameter(const TString& name); //! Fix or release the resonance mass /*! \param [in] parFixed new status of mass */ void fixMass(const Bool_t parFixed) { if (resMass_!=0) { resMass_->fixed(parFixed); } } //! Fix or release the resonance width /*! \param [in] parFixed new status of width */ void fixWidth(const Bool_t parFixed) { if (resWidth_!=0) { resWidth_->fixed(parFixed); } } //! Get the status of resonance mass (fixed or released) /*! \return the status of resonance mass (fixed or released) */ Bool_t fixMass() const { return (resMass_!=0) ? resMass_->fixed() : kTRUE; } //! Get the status of resonance width (fixed or released) /*! \return the status of resonance width (fixed or released) */ Bool_t fixWidth() const { return (resWidth_!=0) ? resWidth_->fixed() : kTRUE; } //! Set the spin formalism to be used /*! \param [in] spinType the spin formalism */ void setSpinType(const LauSpinType spinType) {spinType_ = spinType;} //! Set the form factor model and parameters /*! \param [in] resFactor the barrier factor for the resonance decay \param [in] parFactor the barrier factor for the parent decay */ void setBarrierRadii(LauBlattWeisskopfFactor* resFactor, LauBlattWeisskopfFactor* parFactor) { resBWFactor_ = resFactor; parBWFactor_ = parFactor; } //! Fix or release the Blatt-Weisskopf barrier radii void fixBarrierRadii(const Bool_t fixResRadius, const Bool_t fixParRadius); //! Get the status of resonance barrier radius (fixed or released) Bool_t fixResRadius() const; //! Get the status of parent barrier radius (fixed or released) Bool_t fixParRadius() const; //! Get the radius of the resonance barrier factor Double_t getResRadius() const; //! Get the radius of the parent barrier factor Double_t getParRadius() const; protected: //! Get the name of the parent particle TString getNameParent() const; //! Get the name of the first daughter of the resonance TString getNameDaug1() const; //! Get the name of the second daughter of the resonance TString getNameDaug2() const; //! Get the name of the daughter that does not originate form the resonance TString getNameBachelor() const; //! Get the parent particle mass Double_t getMassParent() const; //! Get the mass of daughter 1 Double_t getMassDaug1() const; //! Get the mass of daughter 2 Double_t getMassDaug2() const; //! Get the mass of the bachelor daughter Double_t getMassBachelor() const; //! Get the Charge of the parent particle Int_t getChargeParent() const; //! Get the charge of daughter 1 Int_t getChargeDaug1() const; //! Get the charge of daughter 2 Int_t getChargeDaug2() const; //! Get the charge of the bachelor daughter Int_t getChargeBachelor() const; //! Get the current value of the daughter momentum in the resonance rest frame Double_t getQ() const {return q_;} //! Get the current value of the bachelor momentum in the resonance rest frame Double_t getP() const {return p_;} //! Get the current value of the bachelor momentum in the parent rest frame Double_t getPstar() const {return pstar_;} //! Get the current value of the full spin-dependent covariant factor Double_t getCovFactor() const {return covFactor_;} //! Get the centrifugal barrier for the parent decay LauBlattWeisskopfFactor* getParBWFactor() {return parBWFactor_;} const LauBlattWeisskopfFactor* getParBWFactor() const {return parBWFactor_;} //! Get the centrifugal barrier for the resonance decay LauBlattWeisskopfFactor* getResBWFactor() {return resBWFactor_;} const LauBlattWeisskopfFactor* getResBWFactor() const {return resBWFactor_;} //! Access the resonance info object LauResonanceInfo* getResInfo() const {return resInfo_;} //! Access the daughters object const LauDaughters* getDaughters() const {return daughters_;} //! Calculate the amplitude spin term using the Zemach tensor formalism /*! \param [in] pProd the momentum factor (either q * p or q * pstar) */ Double_t calcZemachSpinFactor( const Double_t pProd ) const; //! Calculate the amplitude spin term using the covariant tensor formalism /*! \param [in] pProd the momentum factor (q * pstar) */ Double_t calcCovSpinFactor( const Double_t pProd ); //! Calculate the spin-dependent covariant factor /*! \param [in] erm E_ij in the parent rest-frame divided by m_ij (equivalent to sqrt(1 + p^2/mParent^2)) */ void calcCovFactor( const Double_t erm ); //! Calculate the Legendre polynomial for the spin factor /*! Uses the current-event value of cosHel_ */ Double_t calcLegendrePoly() const; //! Calculate the Legendre polynomial for the spin factor (specifying the cosHel value) /*! \param [in] cosHel the cosine of the helicity angle */ Double_t calcLegendrePoly( const Double_t cosHel ); //! Complex resonant amplitude /*! \param [in] mass appropriate invariant mass for the resonance \param [in] spinTerm spin term */ virtual LauComplex resAmp(Double_t mass, Double_t spinTerm) = 0; //! Clear list of floating parameters void clearFloatingParameters() { resParameters_.clear(); } //! Add parameter to the list of floating parameters /*! \param [in] param the parameter to be added to the list */ void addFloatingParameter( LauParameter* param ); //! Access the list of floating parameters std::vector& getParameters() { return resParameters_; } private: //! Copy constructor (not implemented) LauAbsResonance(const LauAbsResonance& rhs); //! Copy assignment operator (not implemented) LauAbsResonance& operator=(const LauAbsResonance& rhs); //! Information on the resonance LauResonanceInfo* resInfo_; //! Information on the particles const LauDaughters* daughters_; //! Parent name TString nameParent_; //! Daughter 1 name TString nameDaug1_; //! Daughter 2 name TString nameDaug2_; //! Bachelor name TString nameBachelor_; //! Parent charge Int_t chargeParent_; //! Daughter 1 charge Int_t chargeDaug1_; //! Daughter 2 charge Int_t chargeDaug2_; //! Bachelor charge Int_t chargeBachelor_; //! Parent mass Double_t massParent_; //! Daughter 1 mass Double_t massDaug1_; //! Daughter 2 mass Double_t massDaug2_; // Bachelor mass Double_t massBachelor_; //! Resonance name TString resName_; //! Resonance name with illegal characters removed TString sanitisedName_; //! Resonance mass LauParameter* resMass_; //! Resonance width LauParameter* resWidth_; //! All parameters of the resonance std::vector resParameters_; //! Resonance spin Int_t resSpin_; //! Resonance charge Int_t resCharge_; //! DP axis identifier Int_t resPairAmpInt_; //! Blatt Weisskopf barrier for parent decay LauBlattWeisskopfFactor* parBWFactor_; //! Blatt Weisskopf barrier for resonance decay LauBlattWeisskopfFactor* resBWFactor_; //! Spin formalism LauSpinType spinType_; //! Boolean to flip helicity Bool_t flipHelicity_; //! Boolean to ignore the momentum factors in both the spin factor and the mass-dependent width Bool_t ignoreMomenta_; //! Boolean to set the spinTerm to unity always Bool_t ignoreSpin_; //! Boolean to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width Bool_t ignoreBarrierScaling_; // Event kinematics information //! Invariant mass Double_t mass_; //! Helicity angle cosine Double_t cosHel_; //! Daughter momentum in resonance rest frame Double_t q_; //! Bachelor momentum in resonance rest frame Double_t p_; //! Bachelor momentum in parent rest frame Double_t pstar_; //! Covariant factor /*! sqrt(1 + z*z), where z = p / mParent Can also be expressed as E_ij in the parent rest-frame divided by m_ij - indeed this is how LauKinematics calculates it. \see LauKinematics::getcov12 \see LauKinematics::getcov13 \see LauKinematics::getcov23 */ Double_t erm_; //! Covariant factor (full spin-dependent expression) Double_t covFactor_; ClassDef(LauAbsResonance,0) // Abstract resonance class }; #endif diff --git a/inc/Laura++_LinkDef.h b/inc/Laura++_LinkDef.h index 4f3296b..0e5441a 100644 --- a/inc/Laura++_LinkDef.h +++ b/inc/Laura++_LinkDef.h @@ -1,155 +1,156 @@ /* Copyright 2013 University of Warwick Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except in compliance with the License. You may obtain a copy of the License at http://www.apache.org/licenses/LICENSE-2.0 Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License. */ /* Laura++ package authors: John Back Paul Harrison Thomas Latham */ #ifdef __CINT__ #pragma link off all globals; #pragma link off all classes; #pragma link off all functions; #pragma link C++ class Lau1DCubicSpline+; #pragma link C++ class Lau1DHistPdf+; #pragma link C++ class Lau2DAbsDP+; #pragma link C++ class Lau2DAbsDPPdf+; #pragma link C++ class Lau2DAbsHistDP+; #pragma link C++ class Lau2DAbsHistDPPdf+; #pragma link C++ class Lau2DCubicSpline+; #pragma link C++ class Lau2DHistDP+; #pragma link C++ class Lau2DHistDPPdf+; #pragma link C++ class Lau2DHistPdf+; #pragma link C++ class Lau2DSplineDP+; #pragma link C++ class Lau2DSplineDPPdf+; #pragma link C++ class LauAbsBkgndDPModel+; #pragma link C++ class LauAbsCoeffSet+; #pragma link C++ class LauAbsEffModel+; #pragma link C++ class LauAbsFitter+; #pragma link C++ class LauAbsFitModel+; #pragma link C++ class LauAbsIncohRes+; #pragma link C++ class LauAbsModIndPartWave+; #pragma link C++ class LauAbsPdf+; #pragma link C++ class LauAbsResonance+; #pragma link C++ class LauAbsRValue+; #pragma link C++ class LauArgusPdf+; #pragma link C++ class LauAsymmCalc+; #pragma link C++ class LauBelleCPCoeffSet+; #pragma link C++ class LauBelleNR+; +#pragma link C++ class LauRescatterThreshold+; #pragma link C++ class LauBelleSymNR+; #pragma link C++ class LauBifurcatedGaussPdf+; #pragma link C++ class LauBkgndDPModel+; #pragma link C++ class LauBlattWeisskopfFactor+; #pragma link C++ class LauBlind+; #pragma link C++ class LauBreitWignerRes+; #pragma link C++ class LauCacheData+; #pragma link C++ class LauCalcChiSq+; #pragma link C++ class LauCartesianCPCoeffSet+; #pragma link C++ class LauCartesianGammaCPCoeffSet+; #pragma link C++ class LauChebychevPdf+; #pragma link C++ class LauCleoCPCoeffSet+; #pragma link C++ class LauComplex+; #pragma link C++ class LauCPFitModel+; #pragma link C++ class LauCruijffPdf+; #pragma link C++ class LauCrystalBallPdf+; #pragma link C++ class LauDabbaRes+; #pragma link C++ class LauDatabasePDG+; #pragma link C++ class LauDaughters+; #pragma link C++ class LauDPDepBifurGaussPdf+; #pragma link C++ class LauDPDepCruijffPdf+; #pragma link C++ class LauDPDepGaussPdf+; #pragma link C++ class LauDPDepMapPdf+; #pragma link C++ class LauDPDepSumPdf+; #pragma link C++ class LauEffModel+; #pragma link C++ class LauEFKLLMRes+; #pragma link C++ class LauEmbeddedData+; #pragma link C++ class LauExponentialPdf+; #pragma link C++ class LauFitDataTree+; #pragma link C++ class LauFitNtuple+; #pragma link C++ class LauFitter+; #pragma link C++ class LauFitObject+; #pragma link C++ class LauFlatteRes+; #pragma link C++ class LauFlatNR+; #pragma link C++ class LauFormulaPar+; #pragma link C++ class LauGaussIncohRes+; #pragma link C++ class LauGaussPdf+; #pragma link C++ class LauGenNtuple+; #pragma link C++ class LauGounarisSakuraiRes+; #pragma link C++ class LauIntegrals+; #pragma link C++ class LauDPPartialIntegralInfo+; #pragma link C++ class LauIsobarDynamics+; #pragma link C++ class LauKappaRes+; #pragma link C++ class LauKinematics+; #pragma link C++ class LauKMatrixProdPole+; #pragma link C++ class LauKMatrixProdSVP+; #pragma link C++ class LauKMatrixPropagator+; #pragma link C++ class LauKMatrixPropFactory+; #pragma link C++ class LauLASSBWRes+; #pragma link C++ class LauLASSNRRes+; #pragma link C++ class LauLASSRes+; #pragma link C++ class LauLinearPdf+; #pragma link C++ class LauMagPhaseCoeffSet+; #pragma link C++ class LauMagPhaseCPCoeffSet+; #pragma link C++ class LauMergeDataFiles+; #pragma link C++ class LauMinuit+; #pragma link C++ class LauModIndPartWaveMagPhase+; #pragma link C++ class LauModIndPartWaveRealImag+; #pragma link C++ class LauNovosibirskPdf+; #pragma link C++ class LauNRAmplitude+; #pragma link C++ class LauParameter+; #pragma link C++ class LauParametricStepFuncPdf+; #pragma link C++ class LauParamFixed+; #pragma link C++ class LauParticlePDG+; #pragma link C++ class LauPolNR+; #pragma link C++ class LauPoleRes+; #pragma link C++ class LauPolarFormFactorNR+; #pragma link C++ class LauPolarFormFactorSymNR+; #pragma link C++ class LauPolarGammaCPCoeffSet+; #pragma link C++ class LauPrint+; #pragma link C++ class LauRealImagCoeffSet+; #pragma link C++ class LauRealImagCPCoeffSet+; #pragma link C++ class LauRealImagGammaCPCoeffSet+; #pragma link C++ class LauRelBreitWignerRes+; #pragma link C++ class LauResonanceInfo+; #pragma link C++ class LauRescatteringRes+; #pragma link C++ class LauRescattering2Res+; #pragma link C++ class LauResonanceMaker+; #pragma link C++ class LauResultsExtractor+; #pragma link C++ class LauRhoOmegaMix+; #ifdef DOLAUROOFITSLAVE #pragma link C++ class LauRooFitSlave+; #endif #pragma link C++ class LauScfMap+; #pragma link C++ class LauSigmaRes+; #pragma link C++ class LauSigmoidPdf+; #pragma link C++ class LauSimpleFitModel+; #pragma link C++ class LauSimFitMaster+; #pragma link C++ class LauSimFitSlave+; #pragma link C++ class LauSPlot+; #pragma link C++ class LauString+; #pragma link C++ class LauSumPdf+; #pragma link C++ class LauTextFileParser+; #pragma link C++ class LauVetoes+; #pragma link C++ class LauWeightedSumEffModel+; #pragma link C++ namespace LauConstants+; #pragma link C++ namespace LauRandom+; #endif diff --git a/src/LauAbsResonance.cc b/src/LauAbsResonance.cc index f956be6..e64116b 100644 --- a/src/LauAbsResonance.cc +++ b/src/LauAbsResonance.cc @@ -1,720 +1,721 @@ /* Copyright 2004 University of Warwick Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except in compliance with the License. You may obtain a copy of the License at http://www.apache.org/licenses/LICENSE-2.0 Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License. */ /* Laura++ package authors: John Back Paul Harrison Thomas Latham */ /*! \file LauAbsResonance.cc \brief File containing implementation of LauAbsResonance class. */ #include #include "TSystem.h" #include "LauAbsResonance.hh" #include "LauConstants.hh" #include "LauDaughters.hh" #include "LauKinematics.hh" #include "LauParameter.hh" #include "LauResonanceInfo.hh" ClassImp(LauAbsResonance) bool LauAbsResonance::isIncoherentModel(LauResonanceModel model) { switch(model) { case BW: case RelBW: case GS: case Flatte: case Sigma: case Kappa: case Dabba: case LASS: case LASS_BW: case LASS_NR: case EFKLLM: case KMatrix: case FlatNR: case NRModel: case BelleNR: + case RescatterThreshold: case PowerLawNR: case BelleSymNR: case BelleSymNRNoInter: case TaylorNR: case PolNR: case Pole: case PolarFFNR: case PolarFFSymNR: case PolarFFSymNRNoInter: case Rescattering: case Rescattering2: case RescatteringNoInter: case MIPW_MagPhase: case MIPW_RealImag: case RhoOmegaMix_GS: case RhoOmegaMix_RBW: case RhoOmegaMix_GS_1: case RhoOmegaMix_RBW_1: break; case GaussIncoh: return true; } return false; } // Constructor LauAbsResonance::LauAbsResonance(LauResonanceInfo* resInfo, const Int_t resPairAmpInt, const LauDaughters* daughters) : resInfo_(resInfo), daughters_(daughters), nameParent_(""), nameDaug1_(""), nameDaug2_(""), nameBachelor_(""), chargeParent_(0), chargeDaug1_(0), chargeDaug2_(0), chargeBachelor_(0), massParent_(0.0), massDaug1_(0.0), massDaug2_(0.0), massBachelor_(0.0), resName_( (resInfo!=0) ? resInfo->getName() : "" ), sanitisedName_( (resInfo!=0) ? resInfo->getSanitisedName() : "" ), resMass_( (resInfo!=0) ? resInfo->getMass() : 0 ), resWidth_( (resInfo!=0) ? resInfo->getWidth() : 0 ), resSpin_( (resInfo!=0) ? resInfo->getSpin() : 0 ), resCharge_( (resInfo!=0) ? resInfo->getCharge() : 0 ), resPairAmpInt_(resPairAmpInt), parBWFactor_(0), resBWFactor_(0), spinType_(Zemach_P), flipHelicity_(kFALSE), ignoreMomenta_(kFALSE), ignoreSpin_(kFALSE), ignoreBarrierScaling_(kFALSE), mass_(0.0), cosHel_(0.0), q_(0.0), p_(0.0), pstar_(0.0), erm_(1.0), covFactor_(1.0) { if ( resInfo == 0 ) { std::cerr << "ERROR in LauAbsResonance constructor : null LauResonanceInfo object provided" << std::endl; gSystem->Exit(EXIT_FAILURE); } if ( daughters_ == 0 ) { std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl; gSystem->Exit(EXIT_FAILURE); } nameParent_ = this->getNameParent(); nameDaug1_ = this->getNameDaug1(); nameDaug2_ = this->getNameDaug2(); nameBachelor_ = this->getNameBachelor(); massParent_ = this->getMassParent(); massDaug1_ = this->getMassDaug1(); massDaug2_ = this->getMassDaug2(); massBachelor_ = this->getMassBachelor(); chargeParent_ = this->getChargeParent(); chargeDaug1_ = this->getChargeDaug1(); chargeDaug2_ = this->getChargeDaug2(); chargeBachelor_ = this->getChargeBachelor(); // check that the total charge adds up to that of the resonance Int_t totalCharge = chargeDaug1_ + chargeDaug2_; if ( (totalCharge != resCharge_) && (resPairAmpInt_ != 0) ) { std::cerr << "ERROR in LauAbsResonance : Total charge of daughters = " << totalCharge << ". Resonance charge = " << resCharge_ << "." << std::endl; gSystem->Exit(EXIT_FAILURE); } } // Constructor LauAbsResonance::LauAbsResonance(const TString& resName, const Int_t resPairAmpInt, const LauDaughters* daughters) : resInfo_(0), daughters_(daughters), nameParent_(""), nameDaug1_(""), nameDaug2_(""), nameBachelor_(""), chargeParent_(0), chargeDaug1_(0), chargeDaug2_(0), chargeBachelor_(0), massParent_(0.0), massDaug1_(0.0), massDaug2_(0.0), massBachelor_(0.0), resName_(resName), sanitisedName_(resName), resMass_(0), resWidth_(0), resSpin_(0), resCharge_(0), resPairAmpInt_(resPairAmpInt), parBWFactor_(0), resBWFactor_(0), spinType_(Zemach_P), flipHelicity_(kFALSE), ignoreMomenta_(kFALSE), ignoreSpin_(kFALSE), ignoreBarrierScaling_(kFALSE), mass_(0.0), cosHel_(0.0), q_(0.0), p_(0.0), pstar_(0.0), erm_(1.0), covFactor_(1.0) { if ( daughters_ == 0 ) { std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl; gSystem->Exit(EXIT_FAILURE); } nameParent_ = this->getNameParent(); nameDaug1_ = this->getNameDaug1(); nameDaug2_ = this->getNameDaug2(); nameBachelor_ = this->getNameBachelor(); massParent_ = this->getMassParent(); massDaug1_ = this->getMassDaug1(); massDaug2_ = this->getMassDaug2(); massBachelor_ = this->getMassBachelor(); chargeParent_ = this->getChargeParent(); chargeDaug1_ = this->getChargeDaug1(); chargeDaug2_ = this->getChargeDaug2(); chargeBachelor_ = this->getChargeBachelor(); // Since we haven't been provided with a LauResonanceInfo object we can just // set the change of the resonance to be the sum of the daughter charges resCharge_ = chargeDaug1_ + chargeDaug2_; } // Destructor LauAbsResonance::~LauAbsResonance() { } LauComplex LauAbsResonance::amplitude(const LauKinematics* kinematics) { // Use LauKinematics interface for amplitude // For resonance made from tracks i, j, we need the momenta // of tracks i and k in the i-j rest frame for spin helicity calculations // in the Zemach tensor formalism. // Also need the momentum of track k in the parent rest-frame for // calculation of the Blatt-Weisskopf factors. mass_ = 0.0; cosHel_ = 0.0; q_ = 0.0; p_ = 0.0; pstar_ = 0.0; erm_ = 1.0; covFactor_ = 1.0; if (resPairAmpInt_ == 1) { mass_ = kinematics->getm23(); cosHel_ = kinematics->getc23(); q_ = kinematics->getp2_23(); p_ = kinematics->getp1_23(); pstar_ = kinematics->getp1_Parent(); erm_ = kinematics->getcov23(); } else if (resPairAmpInt_ == 2) { mass_ = kinematics->getm13(); cosHel_ = kinematics->getc13(); q_ = kinematics->getp1_13(); p_ = kinematics->getp2_13(); pstar_ = kinematics->getp2_Parent(); erm_ = kinematics->getcov13(); } else if (resPairAmpInt_ == 3) { mass_ = kinematics->getm12(); cosHel_ = kinematics->getc12(); q_ = kinematics->getp1_12(); p_ = kinematics->getp3_12(); pstar_ = kinematics->getp3_Parent(); erm_ = kinematics->getcov12(); } else { std::cerr << "ERROR in LauAbsResonance::amplitude : Nonsense setup of resPairAmp array." << std::endl; gSystem->Exit(EXIT_FAILURE); } if (this->flipHelicity()) { cosHel_ *= -1.0; } if (this->ignoreMomenta()) { q_ = 1.0; p_ = 1.0; pstar_ = 1.0; erm_ = 1.0; } // Calculate the spin factors Double_t spinTerm(1.0); Double_t pProd(1.0); if (!this->ignoreSpin()) { switch ( this->getSpinType() ) { case Zemach_P: pProd = q_*p_; spinTerm = this->calcZemachSpinFactor( pProd ); break; case Zemach_Pstar: pProd = q_*pstar_; spinTerm = this->calcZemachSpinFactor( pProd ); break; case Covariant: pProd = q_*pstar_; spinTerm = this->calcCovSpinFactor( pProd ); break; case Legendre: spinTerm = this->calcLegendrePoly(); break; } } // Calculate the full amplitude LauComplex resAmplitude = this->resAmp(mass_, spinTerm); return resAmplitude; } void LauAbsResonance::calcCovFactor( const Double_t erm ) { if (resSpin_ == 0) { covFactor_ = 1.0; } else if (resSpin_ == 1) { covFactor_ = erm; } else if (resSpin_ == 2) { covFactor_ = erm*erm + 0.5; } else if (resSpin_ == 3) { covFactor_ = erm*(erm*erm + 1.5); } else if (resSpin_ == 4) { covFactor_ = (8.*erm*erm*erm*erm + 24.*erm*erm + 3.)/35.; } else if (resSpin_ > 4) { std::cerr << "WARNING in LauAbsResonance::calcCovFactor : covariant spin factor cannot (yet) be fully calculated for spin >= 5" << std::endl; std::cerr << " : the function of sqrt(1 + (p/mParent)^2) part will be missing" << std::endl; covFactor_ = 1.0; } } Double_t LauAbsResonance::calcCovSpinFactor( const Double_t pProd ) { if (resSpin_ == 0) { covFactor_ = 1.0; return 1.0; } // Covariant spin factor is (p* q)^L * f_L(erm) * P_L(cosHel) Double_t spinFactor(pProd); for ( Int_t i(1); i < resSpin_; ++i ) { spinFactor *= pProd; } this->calcCovFactor( erm_ ); spinFactor *= covFactor_; spinFactor *= this->calcLegendrePoly(); return spinFactor; } Double_t LauAbsResonance::calcZemachSpinFactor( const Double_t pProd ) const { // Calculate the spin factors // // These are calculated as follows // // -2^j * (q*p)^j * cj * Pj(cosHel) // // where Pj(coshHel) is the jth order Legendre polynomial and // // cj = j! / (2j-1)!! if (resSpin_ == 0) { return 1.0; } Double_t spinFactor(pProd); for ( Int_t i(1); i < resSpin_; ++i ) { spinFactor *= pProd; } spinFactor *= this->calcLegendrePoly(); return spinFactor; } Double_t LauAbsResonance::calcLegendrePoly( const Double_t cosHel ) { cosHel_ = cosHel; return this->calcLegendrePoly(); } Double_t LauAbsResonance::calcLegendrePoly() const { Double_t legPol = 1.0; if (resSpin_ == 1) { // Calculate vector resonance Legendre polynomial legPol = -2.0*cosHel_; } else if (resSpin_ == 2) { // Calculate tensor resonance Legendre polynomial legPol = 4.0*(3.0*cosHel_*cosHel_ - 1.0)/3.0; } else if (resSpin_ == 3) { // Calculate spin 3 resonance Legendre polynomial legPol = -8.0*(5.0*cosHel_*cosHel_*cosHel_ - 3.0*cosHel_)/5.0; } else if (resSpin_ == 4) { // Calculate spin 4 resonance Legendre polynomial legPol = 16.0*(35.0*cosHel_*cosHel_*cosHel_*cosHel_ - 30.0*cosHel_*cosHel_ + 3.0)/35.0; } else if (resSpin_ == 5) { // Calculate spin 5 resonance Legendre polynomial legPol = -32.0*(63.0*cosHel_*cosHel_*cosHel_*cosHel_*cosHel_ - 70.0*cosHel_*cosHel_*cosHel_ + 15.0*cosHel_)/63.0; } else if (resSpin_ > 5) { std::cerr << "WARNING in LauAbsResonance::calcLegendrePoly : Legendre polynomials not (yet) implemented for spin > 5" << std::endl; } return legPol; } void LauAbsResonance::changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin) { if (newMass > 0.0) { resMass_->valueAndRange(newMass,0.0,3.0*newMass); resMass_->initValue(newMass); resMass_->genValue(newMass); std::cout << "INFO in LauAbsResonance::changeResonance : Setting mass to " << resMass_->value() << std::endl; } if (newWidth > 0.0) { resWidth_->valueAndRange(newWidth,0.0,3.0*newWidth); resWidth_->initValue(newWidth); resWidth_->genValue(newWidth); std::cout << "INFO in LauAbsResonance::changeResonance : Setting width to " << resWidth_->value() << std::endl; } if (newSpin > -1) { resSpin_ = newSpin; std::cout << "INFO in LauAbsResonance::changeResonance : Setting spin to " << resSpin_ << std::endl; } } void LauAbsResonance::changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius) { if ( resRadius >= 0.0 && resBWFactor_ != 0 ) { LauParameter* resBWRadius = resBWFactor_->getRadiusParameter(); resBWRadius->value(resRadius); resBWRadius->initValue(resRadius); resBWRadius->genValue(resRadius); std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting resonance factor radius to " << resBWRadius->value() << std::endl; } if ( parRadius >= 0.0 && parBWFactor_ != 0 ) { LauParameter* parBWRadius = parBWFactor_->getRadiusParameter(); parBWRadius->value(parRadius); parBWRadius->initValue(parRadius); parBWRadius->genValue(parRadius); std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting parent factor radius to " << parBWRadius->value() << std::endl; } } void LauAbsResonance::setResonanceParameter(const TString& name, const Double_t value) { //This function should always be overwritten if needed in classes inheriting from LauAbsResonance. std::cerr << "WARNING in LauAbsResonance::setResonanceParameter : Unable to set parameter \"" << name << "\" to value: " << value << "." << std::endl; } void LauAbsResonance::floatResonanceParameter(const TString& name) { //This function should always be overwritten if needed in classes inheriting from LauAbsResonance. std::cerr << "WARNING in LauAbsResonance::floatResonanceParameter : Unable to release parameter \"" << name << "\"." << std::endl; } LauParameter* LauAbsResonance::getResonanceParameter(const TString& name) { //This function should always be overwritten if needed in classes inheriting from LauAbsResonance. std::cerr << "WARNING in LauAbsResonance::getResonanceParameter : Unable to get parameter \"" << name << "\"." << std::endl; return 0; } void LauAbsResonance::addFloatingParameter( LauParameter* param ) { if ( param == 0 ) { return; } if ( param->clone() ) { resParameters_.push_back( param->parent() ); } else { resParameters_.push_back( param ); } } void LauAbsResonance::fixBarrierRadii(const Bool_t fixResRad, const Bool_t fixParRad) { if ( resBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : resonance barrier factor not present, cannot fix/float it" << std::endl; return; } if ( parBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : parent barrier factor not present, cannot fix/float it" << std::endl; return; } LauParameter* resBWRadius = resBWFactor_->getRadiusParameter(); resBWRadius->fixed(fixResRad); LauParameter* parBWRadius = parBWFactor_->getRadiusParameter(); parBWRadius->fixed(fixParRad); } Bool_t LauAbsResonance::fixResRadius() const { if ( resBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::fixResRadius : resonance barrier factor not present" << std::endl; return kTRUE; } LauParameter* bwRadius = resBWFactor_->getRadiusParameter(); return bwRadius->fixed(); } Bool_t LauAbsResonance::fixParRadius() const { if ( parBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::fixParRadius : parent barrier factor not present" << std::endl; return kTRUE; } LauParameter* bwRadius = parBWFactor_->getRadiusParameter(); return bwRadius->fixed(); } Double_t LauAbsResonance::getResRadius() const { if ( resBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::getResRadius : resonance barrier factor not present" << std::endl; return -1.0; } LauParameter* bwRadius = resBWFactor_->getRadiusParameter(); return bwRadius->unblindValue(); } Double_t LauAbsResonance::getParRadius() const { if ( parBWFactor_ == 0 ) { std::cerr << "WARNING in LauAbsResonance::getParRadius : parent barrier factor not present" << std::endl; return -1.0; } LauParameter* bwRadius = parBWFactor_->getRadiusParameter(); return bwRadius->unblindValue(); } Double_t LauAbsResonance::getMassParent() const { // Get the parent mass Double_t mass(LauConstants::mB); if (daughters_) { mass = daughters_->getMassParent(); } return mass; } Double_t LauAbsResonance::getMassDaug1() const { // Get the daughter mass Double_t mass(LauConstants::mPi); if (daughters_) { if (resPairAmpInt_ == 1) { mass = daughters_->getMassDaug2(); } else if (resPairAmpInt_ == 2) { mass = daughters_->getMassDaug1(); } else if (resPairAmpInt_ == 3) { mass = daughters_->getMassDaug1(); } } return mass; } Double_t LauAbsResonance::getMassDaug2() const { // Get the daughter mass Double_t mass(LauConstants::mPi); if (daughters_) { if (resPairAmpInt_ == 1) { mass = daughters_->getMassDaug3(); } else if (resPairAmpInt_ == 2) { mass = daughters_->getMassDaug3(); } else if (resPairAmpInt_ == 3) { mass = daughters_->getMassDaug2(); } } return mass; } Double_t LauAbsResonance::getMassBachelor() const { // Get the bachelor mass Double_t mass(LauConstants::mPi); if (daughters_) { if (resPairAmpInt_ == 1) { mass = daughters_->getMassDaug1(); } else if (resPairAmpInt_ == 2) { mass = daughters_->getMassDaug2(); } else if (resPairAmpInt_ == 3) { mass = daughters_->getMassDaug3(); } } return mass; } Int_t LauAbsResonance::getChargeParent() const { // Get the parent charge Int_t charge(0); if (daughters_) { charge = daughters_->getChargeParent(); } return charge; } Int_t LauAbsResonance::getChargeDaug1() const { // Get the daughter charge Int_t charge(0); if (daughters_) { if (resPairAmpInt_ == 1) { charge = daughters_->getChargeDaug2(); } else if (resPairAmpInt_ == 2) { charge = daughters_->getChargeDaug1(); } else if (resPairAmpInt_ == 3) { charge = daughters_->getChargeDaug1(); } } return charge; } Int_t LauAbsResonance::getChargeDaug2() const { // Get the daughter charge Int_t charge(0); if (daughters_) { if (resPairAmpInt_ == 1) { charge = daughters_->getChargeDaug3(); } else if (resPairAmpInt_ == 2) { charge = daughters_->getChargeDaug3(); } else if (resPairAmpInt_ == 3) { charge = daughters_->getChargeDaug2(); } } return charge; } Int_t LauAbsResonance::getChargeBachelor() const { // Get the bachelor charge Int_t charge(0); if (daughters_) { if (resPairAmpInt_ == 1) { charge = daughters_->getChargeDaug1(); } else if (resPairAmpInt_ == 2) { charge = daughters_->getChargeDaug2(); } else if (resPairAmpInt_ == 3) { charge = daughters_->getChargeDaug3(); } } return charge; } TString LauAbsResonance::getNameParent() const { // Get the parent name TString name(""); if (daughters_) { name = daughters_->getNameParent(); } return name; } TString LauAbsResonance::getNameDaug1() const { // Get the daughter name TString name(""); if (daughters_) { if (resPairAmpInt_ == 1) { name = daughters_->getNameDaug2(); } else if (resPairAmpInt_ == 2) { name = daughters_->getNameDaug1(); } else if (resPairAmpInt_ == 3) { name = daughters_->getNameDaug1(); } } return name; } TString LauAbsResonance::getNameDaug2() const { // Get the daughter name TString name(""); if (daughters_) { if (resPairAmpInt_ == 1) { name = daughters_->getNameDaug3(); } else if (resPairAmpInt_ == 2) { name = daughters_->getNameDaug3(); } else if (resPairAmpInt_ == 3) { name = daughters_->getNameDaug2(); } } return name; } TString LauAbsResonance::getNameBachelor() const { // Get the bachelor name TString name(""); if (daughters_) { if (resPairAmpInt_ == 1) { name = daughters_->getNameDaug1(); } else if (resPairAmpInt_ == 2) { name = daughters_->getNameDaug2(); } else if (resPairAmpInt_ == 3) { name = daughters_->getNameDaug3(); } } return name; }