Changeset View
Changeset View
Standalone View
Standalone View
src/LauRescatterThreshold.cc
- This file was added.
/* | |||||
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 LauRescatterThreshold.cc | |||||
\brief File containing implementation of LauRescatterThreshold class. | |||||
*/ | |||||
#include <iostream> | |||||
#include "LauRescatterThreshold.hh" | |||||
#include "LauDaughters.hh" | |||||
#include "LauParameter.hh" | |||||
#include "LauResonanceInfo.hh" | |||||
ClassImp(LauRescatterThreshold) | |||||
LauRescatterThreshold::LauRescatterThreshold(LauResonanceInfo* resInfo, const Int_t resPairAmpInt, const LauDaughters* daughters) : | |||||
LauAbsResonance(resInfo, resPairAmpInt, daughters), | |||||
_c11(0), | |||||
_c12(0), | |||||
_rescatterChannel_massDaug1(0), | |||||
_rescatterChannel_massDaug2(0), | |||||
forceLegendre_(kTRUE) | |||||
{ | |||||
TString baseParName = this->getSanitisedName(); | |||||
TString parName = baseParName+"_c11"; | |||||
_c11 = resInfo->getExtraParameter( parName ); | |||||
if ( _c11 == 0 ) { | |||||
_c11 = new LauParameter( parName, 0.0, -1.0, 1.0, kTRUE ); | |||||
_c11->secondStage(kTRUE); | |||||
resInfo->addExtraParameter( _c11 ); | |||||
} | |||||
parName = baseParName+"_c12"; | |||||
_c12 = resInfo->getExtraParameter( parName ); | |||||
if ( _c12 == 0 ) { | |||||
_c12 = new LauParameter( parName, 0.0, -2.0, 10.0, kTRUE ); | |||||
_c12->secondStage(kTRUE); | |||||
resInfo->addExtraParameter( _c12 ); | |||||
} | |||||
parName = baseParName+"_rescatterChannel_massDaug1"; | |||||
_rescatterChannel_massDaug1 = resInfo->getExtraParameter( parName ); | |||||
if ( _rescatterChannel_massDaug1 == 0 ) { | |||||
_rescatterChannel_massDaug1 = new LauParameter( parName, 0.0, 0.0, 10.0, kTRUE ); | |||||
_rescatterChannel_massDaug1->secondStage(kTRUE); | |||||
resInfo->addExtraParameter( _rescatterChannel_massDaug1 ); | |||||
} | |||||
parName = baseParName+"_rescatterChannel_massDaug2"; | |||||
_rescatterChannel_massDaug2 = resInfo->getExtraParameter( parName ); | |||||
if ( _rescatterChannel_massDaug2 == 0 ) { | |||||
_rescatterChannel_massDaug2 = new LauParameter( parName, 0.0, 0.0, 10.0, kTRUE ); | |||||
_rescatterChannel_massDaug2->secondStage(kTRUE); | |||||
resInfo->addExtraParameter( _rescatterChannel_massDaug2 ); | |||||
} | |||||
} | |||||
LauRescatterThreshold::~LauRescatterThreshold() | |||||
{ | |||||
} | |||||
void LauRescatterThreshold::initialise() | |||||
{ | |||||
_massDaug1 = this->getMassDaug1(); | |||||
_massDaug2 = this->getMassDaug2(); | |||||
} | |||||
LauComplex LauRescatterThreshold::resAmp(Double_t mass, Double_t spinTerm) | |||||
{ | |||||
// This function returns the complex dynamical amplitude for a rescattering channel near-threshold, | |||||
// given the invariant mass . | |||||
// need: | |||||
// - coupling strength: near-threshold channel self-coupling [p0] | |||||
// - coupling strength: near-treshold channel -> observation channel [p1] | |||||
// - phasespace factor: near-threshold channel | |||||
// - mass | |||||
// - masses of near-threshold channel particles (will need to be set) [p2,p3] | |||||
// - masses of observation-channel particles (use daughters) | |||||
//std::cout<< "Params - c11: " << _c11->value() << std::endl; | |||||
//std::cout<< "Params - c12: " << _c12->value() << std::endl; | |||||
//std::cout<< "Params - rescatterChannel_massDaug1: " << _rescatterChannel_massDaug1->value() << std::endl; | |||||
//std::cout<< "Params - rescatterChannel_massDaug2: " << _rescatterChannel_massDaug2->value() << std::endl; | |||||
//std::cout<< "ResProps - massResDaug1: " << _massDaug1 << std::endl; | |||||
//std::cout<< "ResProps - massResDaug2: " << _massDaug2 << std::endl; | |||||
//std::cout<< "Calc - mass: " << mass << std::endl; | |||||
//std::cout<< "Calc - spinTerm: " << spinTerm << std::endl; | |||||
LauComplex phspFactor1 = phspFactor(mass,_rescatterChannel_massDaug1->value(),_rescatterChannel_massDaug2->value()); | |||||
LauComplex phspFactor2 = phspFactor(mass,_massDaug1,_massDaug2); | |||||
//std::cout<< "Calc - phspFactor1: " << phspFactor1 << std::endl; | |||||
//std::cout<< "Calc - phspFactor2: " << phspFactor2 << std::endl; | |||||
LauComplex num = LauComplex(-_c12->value(),0) / phspFactor2 ; | |||||
LauComplex im(0,1); | |||||
LauComplex den = LauComplex( _c11->value(),0) - im*phspFactor1 - im*LauComplex(_c12->value()*_c12->value(),0)/phspFactor2 ; | |||||
LauComplex resAmplitude = num / den ; | |||||
//std::cout<< "Calc - num: " << num << std::endl; | |||||
//std::cout<< "Calc - den: " << den << std::endl; | |||||
//std::cout<< "Calc - amp before scale: " << resAmplitude << std::endl; | |||||
resAmplitude.rescale(spinTerm); | |||||
//std::cout<< "Calc - amp after scale: " << resAmplitude << std::endl; | |||||
return resAmplitude; | |||||
} | |||||
LauComplex LauRescatterThreshold::phspFactor(Double_t mass, Double_t m1, Double_t m2) | |||||
{ | |||||
LauComplex phspFactor(1,0); | |||||
Double_t _fact = (mass*mass - (m1+m2)*(m1+m2)) * (mass*mass - (m1-m2)*(m1-m2)); | |||||
if (_fact>0) phspFactor = LauComplex(std::sqrt(_fact),0); | |||||
else phspFactor = LauComplex(0,std::sqrt(-1*_fact)); | |||||
phspFactor.rescale(1./mass); | |||||
return phspFactor; | |||||
} | |||||
const std::vector<LauParameter*>& LauRescatterThreshold::getFloatingParameters() | |||||
{ | |||||
this->clearFloatingParameters(); | |||||
if ( ! this->fix_c11() ) { | |||||
this->addFloatingParameter( _c11 ); | |||||
} | |||||
if ( ! this->fix_c12() ) { | |||||
this->addFloatingParameter( _c12 ); | |||||
} | |||||
if ( ! this->fix_rescatterChannel_massDaug1() ) { | |||||
this->addFloatingParameter( _rescatterChannel_massDaug1 ); | |||||
} | |||||
if ( ! this->fix_rescatterChannel_massDaug2() ) { | |||||
this->addFloatingParameter( _rescatterChannel_massDaug2 ); | |||||
} | |||||
return this->getParameters(); | |||||
} | |||||
void LauRescatterThreshold::setResonanceParameter(const TString& name, const Double_t value) | |||||
{ | |||||
// Set various parameters for the lineshape | |||||
if (name == "c11") { | |||||
this->set_c11(value); | |||||
std::cout << "INFO in LauRescatterThreshold::setResonanceParameter : Setting parameter c11 = " << this->get_c11() << std::endl; | |||||
} else if (name == "c12") { | |||||
this->set_c12(value); | |||||
std::cout << "INFO in LauRescatterThreshold::setResonanceParameter : Setting parameter c12 = " << this->get_c12() << std::endl; | |||||
} else if (name == "rescatterChannel_massDaug1") { | |||||
this->set_rescatterChannel_massDaug1(value); | |||||
std::cout << "INFO in LauRescatterThreshold::setResonanceParameter : Setting parameter rescatterChannel_massDaug1 = " << this->get_rescatterChannel_massDaug1() << std::endl; | |||||
} else if (name == "rescatterChannel_massDaug2") { | |||||
this->set_rescatterChannel_massDaug2(value); | |||||
std::cout << "INFO in LauRescatterThreshold::setResonanceParameter : Setting parameter rescatterChannel_massDaug2 = " << this->get_rescatterChannel_massDaug2() << std::endl; | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::setResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl; | |||||
} | |||||
} | |||||
void LauRescatterThreshold::floatResonanceParameter(const TString& name) | |||||
{ | |||||
if (name == "c11") { | |||||
if ( _c11->fixed() ) { | |||||
_c11->fixed( kFALSE ); | |||||
this->addFloatingParameter( _c11 ); | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl; | |||||
} | |||||
} else if (name == "c12") { | |||||
if ( _c12->fixed() ) { | |||||
_c12->fixed( kFALSE ); | |||||
this->addFloatingParameter( _c12 ); | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl; | |||||
} | |||||
} else if (name == "rescatterChannel_massDaug1") { | |||||
if ( _rescatterChannel_massDaug1->fixed() ) { | |||||
_rescatterChannel_massDaug1->fixed( kFALSE ); | |||||
this->addFloatingParameter( _rescatterChannel_massDaug1 ); | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl; | |||||
} | |||||
} else if (name == "rescatterChannel_massDaug2") { | |||||
if ( _rescatterChannel_massDaug2->fixed() ) { | |||||
_rescatterChannel_massDaug2->fixed( kFALSE ); | |||||
this->addFloatingParameter( _rescatterChannel_massDaug2 ); | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl; | |||||
} | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl; | |||||
} | |||||
} | |||||
LauParameter* LauRescatterThreshold::getResonanceParameter(const TString& name) | |||||
{ | |||||
if (name == "c11") { | |||||
return _c11; | |||||
} else if (name == "c12") { | |||||
return _c12; | |||||
} else if (name == "rescatterChannel_massDaug1") { | |||||
return _rescatterChannel_massDaug1; | |||||
} else if (name == "rescatterChannel_massDaug2") { | |||||
return _rescatterChannel_massDaug2; | |||||
} else { | |||||
std::cerr << "WARNING in LauRescatterThreshold::getResonanceParameter: Parameter name not reconised." << std::endl; | |||||
return 0; | |||||
} | |||||
} | |||||
void LauRescatterThreshold::set_c11(const Double_t c11) | |||||
{ | |||||
_c11->value( c11 ); | |||||
_c11->genValue( c11 ); | |||||
_c11->initValue( c11 ); | |||||
} | |||||
void LauRescatterThreshold::set_c12(const Double_t c12) | |||||
{ | |||||
_c12->value( c12 ); | |||||
_c12->genValue( c12 ); | |||||
_c12->initValue( c12 ); | |||||
} | |||||
void LauRescatterThreshold::set_rescatterChannel_massDaug1(const Double_t rescatterChannel_massDaug1) | |||||
{ | |||||
_rescatterChannel_massDaug1->value( rescatterChannel_massDaug1 ); | |||||
_rescatterChannel_massDaug1->genValue( rescatterChannel_massDaug1 ); | |||||
_rescatterChannel_massDaug1->initValue( rescatterChannel_massDaug1 ); | |||||
} | |||||
void LauRescatterThreshold::set_rescatterChannel_massDaug2(const Double_t rescatterChannel_massDaug2) | |||||
{ | |||||
_rescatterChannel_massDaug2->value( rescatterChannel_massDaug2 ); | |||||
_rescatterChannel_massDaug2->genValue( rescatterChannel_massDaug2 ); | |||||
_rescatterChannel_massDaug2->initValue( rescatterChannel_massDaug2 ); | |||||
} | |||||
No newline at end of file |