Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
//************************************************************************//
// //
// Copyright 2013 Bertram Kopf (bertram@ep1.rub.de) //
// Julian Pychy (julian@ep1.rub.de) //
// - Ruhr-Universität Bochum //
// //
// This file is part of Pawian. //
// //
// Pawian is free software: you can redistribute it and/or modify //
// it under the terms of the GNU General Public License as published by //
// the Free Software Foundation, either version 3 of the License, or //
// (at your option) any later version. //
// //
// Pawian is distributed in the hope that it will be useful, //
// but WITHOUT ANY WARRANTY; without even the implied warranty of //
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
// GNU General Public License for more details. //
// //
// You should have received a copy of the GNU General Public License //
// along with Pawian. If not, see <http://www.gnu.org/licenses/>. //
// //
//************************************************************************//
// KMatrixDynamics class definition file. -*- C++ -*-
// Copyright 2013 Bertram Kopf
#include <getopt.h>
#include <fstream>
#include <string>
#include "PwaUtils/KMatrixDynamics.hh"
#include "PwaUtils/XdecAmpRegistry.hh"
#include "PwaUtils/AbsDecay.hh"
#include "PwaUtils/AbsXdecAmp.hh"
#include "PwaUtils/GlobalEnv.hh"
#include "ErrLogger/ErrLogger.hh"
#include "Particle/Particle.hh"
#include "Particle/ParticleTable.hh"
#include "PwaDynamics/FVector.hh"
#include "ConfigParser/KMatrixParser.hh"
#include "PwaDynamics/KMatrixRel.hh"
#include "PwaDynamics/KMatrixRelBg.hh"
#include "PwaDynamics/KPole.hh"
#include "PwaDynamics/KPoleBarrier.hh"
#include "PwaDynamics/PPole.hh"
#include "PwaDynamics/PPoleBarrier.hh"
#include "PwaDynamics/AbsPhaseSpace.hh"
#include "PwaDynamics/PhaseSpaceIsobar.hh"
#include "FitParams/AbsPawianParameters.hh"
KMatrixDynamics::KMatrixDynamics(std::string& name, std::vector<Particle*>& fsParticles, Particle* mother, std::string& pathToConfigParser) :
AbsDynamics(name, fsParticles, mother)
,_kMatName("")
, _orderKMatBg(-1)
,_withKMatAdler(false)
,_currentAdler0(0.)
,_kMatrixParser(new KMatrixParser(pathToConfigParser))
{
init();
Bertram Kopf
committed
_isLdependent=false;
}
KMatrixDynamics::~KMatrixDynamics()
{
}
complex<double> KMatrixDynamics::eval(EvtData* theData, AbsXdecAmp* grandmaAmp, Spin OrbMom){
complex<double> result(0.,0.);
int evtNo=theData->evtNo;
std::string currentKey="default";
if(0!=grandmaAmp) currentKey=_massKey+grandmaAmp->absDec()->massParKey();
if ( _cacheAmps && !_recalcMap.at(currentKey)){
Bertram Kopf
committed
result=_cachedStringMap.at(evtNo).at(currentKey);
theMutex.lock();
result=_fVecMap.at(currentKey)->evalProjMatrix(theData->DoubleString.at(_dynKey), _projectionIndex);
if ( _cacheAmps){
_cachedStringMap[evtNo][currentKey]=result;
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
// bool readFromCachedMap=false;
// theMutex.lock();
// if( _cacheAmps){
// if(_recalculate){
// bool currentEvtAlreadyCached=false;
// std::map<int, std::map<std::string, bool > >::iterator itAlreadyCached=_alreadyCached.find(evtNo);
// if( itAlreadyCached != _alreadyCached.end()){
// std::map<std::string, bool >::iterator itAlreadyCached2= itAlreadyCached->second.find(currentKey);
// if( itAlreadyCached2 != itAlreadyCached->second.end()){
// currentEvtAlreadyCached=itAlreadyCached2->second;
// }
// else{
// _alreadyCached[evtNo][currentKey]=false;
// }
// }
// else{
// _alreadyCached[evtNo][currentKey]=false;
// }
// if(currentEvtAlreadyCached) readFromCachedMap=true;
// }
// else readFromCachedMap=true;
// }
// if ( readFromCachedMap){
// result=_cachedStringMap.at(evtNo).at(currentKey);
// }
// else{
// result=_fVecMap[currentKey]->evalProjMatrix(theData->DoubleString.at(_dynKey), _projectionIndex);
// if ( _cacheAmps){
// _cachedStringMap[evtNo][currentKey]=result;
// _alreadyCached.at(evtNo).at(currentKey)=true;
// }
// }
// theMutex.unlock();
Bertram Kopf
committed
return result;
}
void KMatrixDynamics::fillDefaultParams(std::shared_ptr<AbsPawianParameters> fitPar){
//beta factor for production
std::map<std::string, std::map<std::string, double> >::iterator it1;
for(it1=_currentbFactorMap.begin(); it1!=_currentbFactorMap.end(); ++it1){
std::string theName=it1->first;
std::map<std::string, double>& bFactors = it1->second;
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::string currentName="b_"+_poleNames.at(i);
// std::cout << "currentName: " << currentName << std::endl;
std::string magName=currentName+"Mag";
fitPar->Add(theName+magName, bFactors.at(magName) , 1.);
// fitPar->SetLimits(theName+magName, 0., bFactors.at(magName)+30.);
std::string phiName=currentName+"Phi";
fitPar->Add(theName+phiName, bFactors.at(phiName) , 0.2);
}
}
//pole positions
std::vector<double >::iterator itPoleVec;
for(unsigned int i=0; i<_currentPoleMasses.size(); ++i){
double valMass=_currentPoleMasses.at(i);
double errMass=0.02;
double minMass=valMass-5.*errMass;
if(minMass<0.) minMass=0.;
double maxMass=valMass+5.*errMass;
fitPar->Add(_poleNames.at(i)+"Mass", valMass, errMass);
fitPar->SetLimits(_poleNames.at(i)+"Mass", minMass, maxMass);
}
//g-factors
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::vector<double> currentgFactorVec=_currentgFactorMap.at(i);
for(unsigned int j=0; j<currentgFactorVec.size(); ++j){
std::string currentName=_poleNames.at(i)+_gFactorNames.at(j)+"gFactor";
double currentError=currentgFactorVec.at(j)/3.;
if (currentError<1.e-5) currentError=0.01;
fitPar->Add(currentName, currentgFactorVec.at(j), currentError);
}
}
//k-matrix bg-terms
if(_orderKMatBg>=0){
for(unsigned int i=0; i<=fabs(_orderKMatBg); ++i){
for(unsigned int j=0; j<_phpVecs.size(); ++j){
for(unsigned int k=j; k<_phpVecs.size(); ++k){
std::string currentName=_bgTermNames.at(i).at(j).at(k);
fitPar->Add(currentName, _currentBgTerms.at(i).at(j).at(k), fabs(_currentBgTerms.at(i).at(j).at(k))+0.3);
}
}
}
//Adler-term
if(_withKMatAdler){
fitPar->Add("s0"+_kMatName, _currentAdler0, fabs(_currentAdler0)+0.2);
}
}
}
Bertram Kopf
committed
void KMatrixDynamics::fillParamNameList(){
_paramNameList.clear();
Bertram Kopf
committed
Bertram Kopf
committed
//beta factor for production
std::map<std::string, std::map<std::string, double> >::iterator it1;
for(it1=_currentbFactorMap.begin(); it1!=_currentbFactorMap.end(); ++it1){
Bertram Kopf
committed
std::string theName=it1->first;
std::vector<std::string> currentNameList;
Bertram Kopf
committed
std::map<std::string, double>& bFactors = it1->second;
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::string currentName="b_"+_poleNames.at(i);
// std::cout << "currentName: " << currentName << std::endl;
std::string magName=currentName+"Mag";
_paramNameList.push_back(theName+magName);
currentNameList.push_back(theName+magName);
Bertram Kopf
committed
std::string phiName=currentName+"Phi";
_paramNameList.push_back(theName+phiName);
currentNameList.push_back(theName+phiName);
_paramNameListMap[theName]=currentNameList;
std::map<std::string, std::vector<std::string> >::iterator itNameMap;
//pole positions
std::vector<double >::iterator itPoleVec;
for(unsigned int i=0; i<_currentPoleMasses.size(); ++i){
Bertram Kopf
committed
_paramNameList.push_back(_poleNames.at(i)+"Mass");
for(itNameMap=_paramNameListMap.begin(); itNameMap!=_paramNameListMap.end(); ++itNameMap){
itNameMap->second.push_back(_poleNames.at(i)+"Mass");
}
Bertram Kopf
committed
//g-factors
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::vector<double> currentgFactorVec=_currentgFactorMap.at(i);
for(unsigned int j=0; j<currentgFactorVec.size(); ++j){
std::string currentName=_poleNames.at(i)+_gFactorNames.at(j)+"gFactor";
Bertram Kopf
committed
_paramNameList.push_back(currentName);
for(itNameMap=_paramNameListMap.begin(); itNameMap!=_paramNameListMap.end(); ++itNameMap){
itNameMap->second.push_back(currentName);
}
}
}
Bertram Kopf
committed
//k-matrix bg-terms
if(_orderKMatBg>=0){
for(unsigned int i=0; i<=fabs(_orderKMatBg); ++i){
for(unsigned int j=0; j<_phpVecs.size(); ++j){
Bertram Kopf
committed
for(unsigned int k=j; k<_phpVecs.size(); ++k){
std::string currentName=_bgTermNames.at(i).at(j).at(k);
_paramNameList.push_back(currentName);
for(itNameMap=_paramNameListMap.begin(); itNameMap!=_paramNameListMap.end(); ++itNameMap){
itNameMap->second.push_back(currentName);
}
Bertram Kopf
committed
}
}
}
//Adler-term
if(_withKMatAdler){
Bertram Kopf
committed
_paramNameList.push_back("s0"+_kMatName);
for(itNameMap=_paramNameListMap.begin(); itNameMap!=_paramNameListMap.end(); ++itNameMap){
itNameMap->second.push_back("s0"+_kMatName);
}
}
}
}
bool KMatrixDynamics::checkRecalculation(std::shared_ptr<AbsPawianParameters> fitParNew, std::shared_ptr<AbsPawianParameters> fitParOld){
std::map<std::string, std::vector<std::string> >::iterator itMap;
for(itMap= _paramNameListMap.begin(); itMap !=_paramNameListMap.end(); ++itMap){
_recalcMap.at(itMap->first)=false;
std::vector<std::string>::iterator itStr;
for (itStr=itMap->second.begin(); itStr!=itMap->second.end(); ++itStr){
std::string currentParamName=*itStr;
if(!CheckDoubleEquality(fitParNew->Value(currentParamName), fitParOld->Value(currentParamName))){
_recalcMap.at(itMap->first)=true;
continue;
}
return AbsParamHandler::checkRecalculation(fitParNew, fitParOld);
Bertram Kopf
committed
void KMatrixDynamics::updateFitParams(std::shared_ptr<AbsPawianParameters> fitPar){
//beta factor for production
std::map<std::string, std::map<std::string, double> >::iterator it1;
for(it1=_currentbFactorMap.begin(); it1!=_currentbFactorMap.end(); ++it1){
Bertram Kopf
committed
std::string theName=it1->first;
std::map<std::string, double>::iterator it2;
std::map<std::string, double>& bFactors=it1->second;
for(it2=bFactors.begin(); it2!=bFactors.end(); ++it2){
Bertram Kopf
committed
it2->second=fitPar->Value(it1->first+it2->first);
Bertram Kopf
committed
std::shared_ptr<PVectorRel> currentPVec=_pVecMap.at(it1->first);
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::string currentName="b_"+_poleNames.at(i);
complex<double> currentbFactor=fabs(bFactors.at(currentName+"Mag"))*complex<double>(cos(bFactors.at(currentName+"Phi")), sin(bFactors.at(currentName+"Phi")));
currentPVec->updateBeta(i, currentbFactor);
}
}
//pole positions
std::vector<double >::iterator itPoleVec;
for(unsigned int i=0; i<_currentPoleMasses.size(); ++i){
Bertram Kopf
committed
std::string currentPoleName=_poleNames.at(i)+"Mass";
double currentPoleMass=fitPar->Value(currentPoleName);
_currentPoleMasses.at(i)=currentPoleMass;
_kPoles.at(i)->updatePoleMass(currentPoleMass);
std::map<std::string, std::shared_ptr<PVectorRel> >::iterator itPVec;
for(itPVec=_pVecMap.begin(); itPVec!=_pVecMap.end(); ++itPVec){
itPVec->second->updatePoleMass(i, currentPoleMass);
}
}
//g-factors
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::vector<double>& currentgFactorVec=_currentgFactorMap.at(i);
for(unsigned int j=0; j<currentgFactorVec.size(); ++j){
std::string currentName=_poleNames.at(i)+_gFactorNames.at(j)+"gFactor";
Bertram Kopf
committed
currentgFactorVec.at(j)=fitPar->Value(currentName);
}
_kPoles.at(i)->updategFactors(currentgFactorVec);
std::map<std::string, std::shared_ptr<PVectorRel> >::iterator itPVec;
for(itPVec=_pVecMap.begin(); itPVec!=_pVecMap.end(); ++itPVec){
itPVec->second->updategFactors(i, currentgFactorVec);
}
}
//k-matrix bg-terms
if(_orderKMatBg>=0){
for(unsigned int i=0; i<=fabs(_orderKMatBg); ++i){
for(unsigned int j=0; j<_phpVecs.size(); ++j){
for(unsigned int k=j; k<_phpVecs.size(); ++k){
Bertram Kopf
committed
std::string currentName=_bgTermNames.at(i).at(j).at(k);
double newVal=fitPar->Value(_bgTermNames.at(i).at(j).at(k));
_currentBgTerms.at(i).at(j).at(k)=newVal;
_kMatr->updateBgTerms(i,j,k,newVal);
}
}
}
//Adler-term
if(_withKMatAdler){
Bertram Kopf
committed
_currentAdler0=fitPar->Value("s0"+_kMatName);
_kMatr->updates0Adler(_currentAdler0);
}
}
Bertram Kopf
committed
Bertram Kopf
committed
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
// void KMatrixDynamics::updateFitParams(fitParCol& theParamVal){
// //beta factor for production
// std::map<std::string, std::map<std::string, double> >::iterator it1;
// for(it1=_currentbFactorMap.begin(); it1!=_currentbFactorMap.end(); ++it1){
// std::map<std::string, double>::iterator it2;
// std::map<std::string, double>& bFactors=it1->second;
// for(it2=bFactors.begin(); it2!=bFactors.end(); ++it2){
// it2->second=theParamVal.otherParams.at(it1->first+it2->first);
// }
// std::shared_ptr<PVectorRel> currentPVec=_pVecMap.at(it1->first);
// for(unsigned int i=0; i<_poleNames.size(); ++i){
// std::string currentName="b_"+_poleNames.at(i);
// complex<double> currentbFactor=bFactors.at(currentName+"Mag")*complex<double>(cos(bFactors.at(currentName+"Phi")), sin(bFactors.at(currentName+"Phi")));
// currentPVec->updateBeta(i, currentbFactor);
// }
// }
// //pole positions
// std::vector<double >::iterator itPoleVec;
// for(unsigned int i=0; i<_currentPoleMasses.size(); ++i){
// double currentPoleMass=theParamVal.Masses[_poleNames.at(i)];
// _currentPoleMasses.at(i)=currentPoleMass;
// _kPoles.at(i)->updatePoleMass(currentPoleMass);
// std::map<std::string, std::shared_ptr<PVectorRel> >::iterator itPVec;
// for(itPVec=_pVecMap.begin(); itPVec!=_pVecMap.end(); ++itPVec){
// itPVec->second->updatePoleMass(i, currentPoleMass);
// }
// }
// //g-factors
// for(unsigned int i=0; i<_poleNames.size(); ++i){
// std::vector<double>& currentgFactorVec=_currentgFactorMap.at(i);
// for(unsigned int j=0; j<currentgFactorVec.size(); ++j){
// std::string currentName=_poleNames.at(i)+_gFactorNames.at(j)+"gFactor";
// currentgFactorVec.at(j)=theParamVal.gFactors.at(currentName);
// }
// _kPoles.at(i)->updategFactors(currentgFactorVec);
// std::map<std::string, std::shared_ptr<PVectorRel> >::iterator itPVec;
// for(itPVec=_pVecMap.begin(); itPVec!=_pVecMap.end(); ++itPVec){
// itPVec->second->updategFactors(i, currentgFactorVec);
// }
// }
// //k-matrix bg-terms
// if(_orderKMatBg>=0){
// for(unsigned int i=0; i<=fabs(_orderKMatBg); ++i){
// for(unsigned int j=0; j<_phpVecs.size(); ++j){
// for(unsigned int k=j; k<_phpVecs.size(); ++k){
// double newVal=theParamVal.otherParams.at(_bgTermNames.at(i).at(j).at(k));
// _currentBgTerms.at(i).at(j).at(k)=newVal;
// _kMatr->updateBgTerms(i,j,k,newVal);
// }
// }
// }
// //Adler-term
// if(_withKMatAdler){
// _currentAdler0=theParamVal.otherParams.at("s0"+_kMatName);
// _kMatr->updates0Adler(_currentAdler0);
// }
// }
// }
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
void KMatrixDynamics::addGrandMa(std::shared_ptr<AbsDecay> theDec){
if(0==theDec){
Alert << "Can not add AbsXdecAmp; 0 pointer!!!" << endmsg;
exit(1);
}
std::string theName=_massKey+theDec->massParKey();
std::cout << "addGrandMa:\t" << theName << std::endl;
std::map<std::string, std::shared_ptr<FVector> >::iterator it = _fVecMap.find(theName);
if (it != _fVecMap.end()) return;
std::shared_ptr<PVectorRel> currentPVector=makeNewPVec();
_pVecMap[theName]=currentPVector;
std::vector< std::string>::iterator poleNameIt;
for (poleNameIt=_poleNames.begin(); poleNameIt!=_poleNames.end(); ++poleNameIt){
std::string currentName="b_"+(*poleNameIt);
_currentbFactorMap[theName][currentName+"Mag"]=1.;
_currentbFactorMap[theName][currentName+"Phi"]=0.;
}
std::map<std::string, double>& bFactors = _currentbFactorMap[theName];
for(unsigned int i=0; i<_poleNames.size(); ++i){
std::string currentName="b_"+_poleNames[i];
std::cout << "currentName: " << currentName << std::endl;
complex<double> currentbFactor=bFactors[currentName+"Mag"]*complex<double>(cos(bFactors[currentName+"Phi"]), sin(bFactors[currentName+"Phi"]));
currentPVector->updateBeta(i, currentbFactor);
}
std::shared_ptr<FVector> currentFVector=std::shared_ptr<FVector>(new FVector(_kMatr, currentPVector));
_fVecMap[theName]=currentFVector;
_recalcMap[theName]=true;
}
const std::string& KMatrixDynamics::grandMaKey(AbsXdecAmp* grandmaAmp){
if(0==grandmaAmp) return _grandmaKey;
return grandmaAmp->absDec()->massParKey();
}
void KMatrixDynamics::init(){
_kMatName=_kMatrixParser->keyName();
_orderKMatBg=_kMatrixParser->orderBg();
_withKMatAdler=_kMatrixParser->useAdler();
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
std::vector<std::string> poleNameAndMassVecs=_kMatrixParser->poles();
std::vector<std::string>::iterator itString;
for (itString=poleNameAndMassVecs.begin(); itString!=poleNameAndMassVecs.end(); ++itString){
std::istringstream poleIString(*itString);
std::string currentPoleName;
std::string currentPoleMassStr;
poleIString >> currentPoleName >> currentPoleMassStr;
std::istringstream currentPoleMassiStr(currentPoleMassStr);
double currentValue;
if(!(currentPoleMassiStr >> currentValue)){
Alert << "cannot convert " << currentPoleMassStr << " to a double value" << endmsg;
exit(0);
}
_currentPoleMasses.push_back(currentValue);
_poleNames.push_back(currentPoleName);
}
if(_currentPoleMasses.size()!= _kMatrixParser->noOfPoles()){
Alert << "number of poles != number of pole masses" << endmsg;
exit(0);
}
const std::vector<std::string> gFacStringVec=_kMatrixParser->gFactors();
std::vector<std::string>::const_iterator itStrConst;
for(itStrConst=gFacStringVec.begin(); itStrConst!=gFacStringVec.end(); ++itStrConst){
std::istringstream particles(*itStrConst);
std::string firstParticleName;
std::string secondParticleName;
particles >> firstParticleName >> secondParticleName;
Particle* firstParticle = GlobalEnv::instance()->particleTable()->particle(firstParticleName);
Particle* secondParticle = GlobalEnv::instance()->particleTable()->particle(secondParticleName);
if(0==firstParticle || 0==secondParticle){
Alert << "particle with name: " << firstParticleName <<" or " << secondParticleName << " doesn't exist in pdg-table" << endmsg;
exit(0);
}
std::shared_ptr<AbsPhaseSpace> currentPhp(new PhaseSpaceIsobar(firstParticle->mass(), secondParticle->mass()));
_phpVecs.push_back(currentPhp);
std::string gFactorKey=firstParticleName+secondParticleName;
_gFactorNames.push_back(gFactorKey);
for (int i=0; i<int(_kMatrixParser->noOfPoles()); ++i){
std::string currentPoleName=_poleNames[i];
std::string currentgValueStr;
if(!(particles >> currentgValueStr)){
Alert << "g-factors for pole " << currentPoleName << " does not exist!" << endmsg;
exit(0);
}
std::istringstream currentgValueiStr(currentgValueStr);
double currentGValue;
if (!(currentgValueiStr >> currentGValue)){
Alert << "cannot convert " << currentgValueStr << " to a double value" << endmsg;
exit(0);
}
std::string gFactorKey=firstParticleName+secondParticleName;
_currentgFactorMap[i].push_back(currentGValue);
}
}
std::map<int, std::vector<double> >::iterator itgFac;
for (itgFac=_currentgFactorMap.begin(); itgFac!=_currentgFactorMap.end(); ++itgFac){
std::vector<double> currentgVector=itgFac->second;
std::shared_ptr<KPole> currentPole;
if (_kMatrixParser->useBarrierFactors()) currentPole=std::shared_ptr<KPole>(new KPoleBarrier(currentgVector, _currentPoleMasses.at(itgFac->first), _phpVecs, _kMatrixParser->orbitalMom(), _kMatrixParser->useTruncatedBarrierFactors()));
else currentPole=std::shared_ptr<KPole>(new KPole(currentgVector, _currentPoleMasses.at(itgFac->first)));
_kPoles.push_back(currentPole);
}
if(_orderKMatBg<0) _kMatr=std::shared_ptr<KMatrixRel>(new KMatrixRel(_kPoles,_phpVecs ));
else {
_currentBgTerms.resize(_orderKMatBg+1);
_bgTermNames.resize(_orderKMatBg+1);
for(unsigned int i=0; i<= fabs(_orderKMatBg); ++i){
_currentBgTerms.at(i).resize(_phpVecs.size());
_bgTermNames.at(i).resize(_phpVecs.size());
for(unsigned int j=0; j<_phpVecs.size(); ++j){
_currentBgTerms.at(i).at(j).resize(_phpVecs.size());
_bgTermNames.at(i).at(j).resize(_phpVecs.size());
for(unsigned int k=j; k<_phpVecs.size(); ++k){
_currentBgTerms.at(i).at(j).at(k)=0.;
std::stringstream keyOrderStrStr;
keyOrderStrStr << i << j << k;
std::string keyOrder=keyOrderStrStr.str();
std::string currentName="bg"+keyOrder+_kMatName;
_bgTermNames.at(i).at(j).at(k)=currentName;
}
}
}
_kMatr=std::shared_ptr<KMatrixRel>(new KMatrixRelBg(_kPoles,_phpVecs, _orderKMatBg, _withKMatAdler ));
if(_withKMatAdler){
_currentAdler0=_kMatrixParser->s0Adler();
_kMatr->updates0Adler(_currentAdler0);
_kMatr->updatesnormAdler(_kMatrixParser->snormAdler());
}
}
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
const std::string porjectionParticleNames=_kMatrixParser->projection();
std::istringstream projParticles(porjectionParticleNames);
std::string firstProjParticleName;
std::string secondProjParticleName;
projParticles >> firstProjParticleName >> secondProjParticleName;
std::string projKey=firstProjParticleName+secondProjParticleName;
bool found=false;
for(unsigned int i=0; i<_gFactorNames.size();++i){
if(projKey==_gFactorNames[i]){
_projectionIndex=i;
found=true;
}
}
if (!found){
Alert << "projection index for key " << projKey << " not found" << endmsg;
exit(0);
}
}
std::shared_ptr<PVectorRel> KMatrixDynamics::makeNewPVec(){
vector<std::shared_ptr<PPole> > thePpoles;
complex<double> defaultBeta(1.,0.);
vector<std::shared_ptr<KPole> >::iterator it;
for (it=_kPoles.begin(); it!=_kPoles.end(); ++it){
std::vector<double> currentGFactors=(*it)->gFactors();
std::shared_ptr<PPole> currentPPole;
if (_kMatrixParser->useBarrierFactors()) currentPPole=std::shared_ptr<PPole>(new PPoleBarrier(defaultBeta, currentGFactors, (*it)->poleMass(), _phpVecs, _kMatrixParser->orbitalMom()));
else currentPPole=std::shared_ptr<PPole>(new PPole(defaultBeta, currentGFactors, (*it)->poleMass()));
thePpoles.push_back(currentPPole);
}
std::shared_ptr<PVectorRel> thePVector(new PVectorRel(thePpoles, _phpVecs));
return thePVector;
}