JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELALinearInterpFunc Class Reference

#include <MELALinearInterpFunc.h>

Inheritance diagram for MELALinearInterpFunc:
Collaboration diagram for MELALinearInterpFunc:

Public Types

enum  VerbosityLevel { kSilent, kError, kVerbose }
 
typedef Float_t T
 

Public Member Functions

 MELALinearInterpFunc ()
 
 MELALinearInterpFunc (const char *name, const char *title)
 
 MELALinearInterpFunc (const char *name, const char *title, RooAbsReal &inXVar, const std::vector< T > &inXList, const RooArgList &inFcnList, Bool_t inUseFloor=true, T inFloorEval=1e-15, T inFloorInt=1e-10)
 
 MELALinearInterpFunc (const MELALinearInterpFunc &other, const char *name=0)
 
virtual TObject * clone (const char *newname) const
 
virtual ~MELALinearInterpFunc ()
 
void setVerbosity (VerbosityLevel flag)
 
void setEvalFloor (T val)
 
void setIntFloor (T val)
 
void doFloor (Bool_t flag)
 
void setRangeValidity (const T valmin, const T valmax)
 
Int_t getAnalyticalIntegral (RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
 
Double_t analyticalIntegral (Int_t code, const char *rangeName=0) const
 

Protected Member Functions

unsigned int npoints () const
 
Int_t getWhichBin (const T &val) const
 
T getKappa (const Int_t &bin) const
 
T getTVar (const T &val) const
 
Bool_t testRangeValidity (const T &val) const
 
void cropValueForRange (T &val) const
 
T interpolateFcn (Int_t code, const char *rangeName=0) const
 
Double_t evaluate () const
 

Protected Attributes

VerbosityLevel verbosity
 
Bool_t useFloor
 
T floorEval
 
T floorInt
 
T rangeXmin
 
T rangeXmax
 
std::vector< TXList
 
RooRealProxy theXVar
 
RooListProxy FcnList
 
RooListProxy leafDepsList
 

Detailed Description

Definition at line 13 of file MELALinearInterpFunc.h.

Member Typedef Documentation

◆ T

typedef Float_t MELALinearInterpFunc::T

Definition at line 15 of file MELALinearInterpFunc.h.

Member Enumeration Documentation

◆ VerbosityLevel

Enumerator
kSilent 
kError 
kVerbose 

Definition at line 17 of file MELALinearInterpFunc.h.

17  {
18  kSilent,
19  kError,
20  kVerbose
21  };

Constructor & Destructor Documentation

◆ MELALinearInterpFunc() [1/4]

MELALinearInterpFunc::MELALinearInterpFunc ( )

◆ MELALinearInterpFunc() [2/4]

MELALinearInterpFunc::MELALinearInterpFunc ( const char *  name,
const char *  title 
)

Definition at line 24 of file MELALinearInterpFunc.cc.

27  :
28  RooAbsReal(name, title),
30  useFloor(true), floorEval(0), floorInt(0),
31  rangeXmin(1), rangeXmax(-1),
32  theXVar("theXVar", "theXVar", this),
33  FcnList("FcnList", "FcnList", this),
34  leafDepsList("leafDepsList", "leafDepsList", this)
35 {}

◆ MELALinearInterpFunc() [3/4]

MELALinearInterpFunc::MELALinearInterpFunc ( const char *  name,
const char *  title,
RooAbsReal &  inXVar,
const std::vector< T > &  inXList,
const RooArgList &  inFcnList,
Bool_t  inUseFloor = true,
T  inFloorEval = 1e-15,
T  inFloorInt = 1e-10 
)

◆ MELALinearInterpFunc() [4/4]

MELALinearInterpFunc::MELALinearInterpFunc ( const MELALinearInterpFunc other,
const char *  name = 0 
)

Definition at line 90 of file MELALinearInterpFunc.cc.

93  :
94  RooAbsReal(other, name),
95  verbosity(other.verbosity),
96  useFloor(other.useFloor), floorEval(other.floorEval), floorInt(other.floorInt),
97  rangeXmin(other.rangeXmin), rangeXmax(other.rangeXmax),
98  XList(other.XList),
99  theXVar("theXVar", this, other.theXVar),
100  FcnList("FcnList", this, other.FcnList),
101  leafDepsList("leafDepsList", this, other.leafDepsList)
102 {}

◆ ~MELALinearInterpFunc()

virtual MELALinearInterpFunc::~MELALinearInterpFunc ( )
inlinevirtual

Definition at line 41 of file MELALinearInterpFunc.h.

41 {}

Member Function Documentation

◆ analyticalIntegral()

Double_t MELALinearInterpFunc::analyticalIntegral ( Int_t  code,
const char *  rangeName = 0 
) const

Definition at line 290 of file MELALinearInterpFunc.cc.

290  {
291  Double_t value = interpolateFcn(code, rangeName);
292  if (useFloor && value<floorInt){
293  if (verbosity>=MELALinearInterpFunc::kError) coutE(Integration) << "MELALinearInterpFunc ERROR::MELALinearInterpFunc(" << GetName() << ") integration returned " << value << " for code = " << code << endl;
294  value = floorInt;
295  }
296  if (verbosity==MELALinearInterpFunc::kVerbose){ cout << "MELALinearInterpFunc(" << GetName() << ")::analyticalIntegral = " << value << " for code = " << code << endl; }
297  return value;
298 }

◆ clone()

virtual TObject* MELALinearInterpFunc::clone ( const char *  newname) const
inlinevirtual

Definition at line 40 of file MELALinearInterpFunc.h.

40 { return new MELALinearInterpFunc(*this, newname); }

◆ cropValueForRange()

void MELALinearInterpFunc::cropValueForRange ( T val) const
protected

Definition at line 313 of file MELALinearInterpFunc.cc.

313  {
314  if (testRangeValidity(val)) return;
315  const T* range[2];
316  range[0] = &rangeXmin;
317  range[1] = &rangeXmax;
318  if (val<*(range[0])) val = *(range[0]);
319  if (val>*(range[1])) val = *(range[1]);
320 }

◆ doFloor()

void MELALinearInterpFunc::doFloor ( Bool_t  flag)

Definition at line 107 of file MELALinearInterpFunc.cc.

107 { useFloor = flag; }

◆ evaluate()

Double_t MELALinearInterpFunc::evaluate ( ) const
protected

Definition at line 222 of file MELALinearInterpFunc.cc.

222  {
223  Double_t value = interpolateFcn(0);
224  if (useFloor && value<floorEval){
225  if (verbosity>=MELALinearInterpFunc::kError) coutE(Eval) << "MELALinearInterpFunc ERROR::MELALinearInterpFunc(" << GetName() << ") evaluation returned " << value << " at x = " << theXVar << endl;
226  value = floorEval;
227  }
228  if (verbosity==MELALinearInterpFunc::kVerbose){ cout << "MELALinearInterpFunc(" << GetName() << ")::evaluate = " << value << " at x = " << theXVar << endl; }
229  return value;
230 }

◆ getAnalyticalIntegral()

Int_t MELALinearInterpFunc::getAnalyticalIntegral ( RooArgSet &  allVars,
RooArgSet &  analVars,
const char *  rangeName = 0 
) const

Definition at line 231 of file MELALinearInterpFunc.cc.

231  {
232  if (_forceNumInt) return 0;
233 
234  const Int_t xprime = 101;
235  Int_t code=1;
236 
237  RooArgSet Xdeps, Fcndeps;
238  RooRealVar* rrv_x = dynamic_cast<RooRealVar*>(theXVar.absArg());
239  if (rrv_x==0) theXVar.absArg()->leafNodeServerList(&Xdeps, 0, true);
240 
241  for (int ix=0; ix<FcnList.getSize(); ix++){
242  RooArgSet tmpdeps;
243  RooAbsReal const* fcn = dynamic_cast<const RooAbsReal*>(FcnList.at(ix));
244  fcn->leafNodeServerList(&tmpdeps, 0, true);
245  Fcndeps.add(tmpdeps);
246  if (ix==0){
247  RooArgSet mirrorvars;
248  TIterator* iter = allVars.createIterator();
249  RooAbsArg* var;
250  while ((var = (RooAbsArg*)iter->Next())){
251  RooRealVar* rrv_tmp = dynamic_cast<RooRealVar*>(var);
252  if (rrv_tmp==0) continue;
253  bool isFound=false;
254  if (rrv_x!=0){
255  isFound = (rrv_tmp==rrv_x);
256  }
257  else{
258  TIterator* iter_deps = Xdeps.createIterator();
259  RooAbsArg* var_deps;
260  while ((var_deps = (RooAbsArg*)iter_deps->Next())){
261  RooRealVar* depvar = dynamic_cast<RooRealVar*>(var_deps);
262  if (depvar==0) continue;
263  if (depvar==rrv_tmp){ isFound=true; break; }
264  }
265  delete iter_deps;
266  }
267  if (!isFound) mirrorvars.add(*var);
268  }
269  delete iter;
270 
271  Int_t codetmp = fcn->getAnalyticalIntegral(mirrorvars, analVars, rangeName);
272  if (codetmp!=0) code *= codetmp;
273  }
274  }
275 
276  if (rrv_x!=0){
277  if (Fcndeps.find(*rrv_x)==0){
278  if (matchArgs(allVars, analVars, theXVar)) code *= xprime;
279  }
280  }
281 
282  if (code==1) code=0;
284  cout << "MELALinearInterpFunc(" << GetName() << ")::getAnalyticalIntegral code = " << code << endl;
285  allVars.Print("v");
286  analVars.Print("v");
287  }
288  return code;
289 }

◆ getKappa()

MELALinearInterpFunc::T MELALinearInterpFunc::getKappa ( const Int_t &  bin) const
protected

Definition at line 130 of file MELALinearInterpFunc.cc.

130  {
131  if (bin>=0){
132  vector<MELALinearInterpFunc::T> const* coord=&XList;
133  MELALinearInterpFunc::T diff = 1;
134  if (Int_t(coord->size())>(bin+1)) diff = coord->at(bin+1)-coord->at(bin);
135  if (fabs(diff)>MELALinearInterpFunc::T(0)) return MELALinearInterpFunc::T(1)/diff;
136  else return 0;
137  }
138  else return 1;
139 }

◆ getTVar()

MELALinearInterpFunc::T MELALinearInterpFunc::getTVar ( const T val) const
protected

Definition at line 140 of file MELALinearInterpFunc.cc.

140  {
141  const Int_t bin = getWhichBin(val);
142  vector<MELALinearInterpFunc::T> const* coord=&XList;
144  return (val - coord->at(bin))*kappa;
145 }

◆ getWhichBin()

Int_t MELALinearInterpFunc::getWhichBin ( const T val) const
protected

Definition at line 109 of file MELALinearInterpFunc.cc.

109  {
110  Int_t bin=-1;
111  MELALinearInterpFunc::T valj, valjpo;
112  Int_t np;
113  vector<MELALinearInterpFunc::T> const* coord=&XList;
114  np=npoints();
115 
116  if (np<=1) bin=0;
117  else{
118  valjpo = coord->at(0);
119  for (Int_t j=0; j<np-1; j++){
120  valj = coord->at(j);
121  valjpo = coord->at(j+1);
122  if (val<valjpo && val>=valj){ bin=j; break; }
123  }
124  if (bin==-1 && val>=valjpo) bin=np-2;
125  else if (bin==-1) bin=0;
126  }
127 
128  return bin;
129 }

◆ interpolateFcn()

MELALinearInterpFunc::T MELALinearInterpFunc::interpolateFcn ( Int_t  code,
const char *  rangeName = 0 
) const
protected

Definition at line 147 of file MELALinearInterpFunc.cc.

147  {
148  const Int_t xprime = 101;
149  const bool intAlongX=(code>0 && code%xprime==0);
150 
152  Int_t coderem = (intAlongX ? code/xprime : code);
153  if (coderem==1) coderem=0;
154 
156  if (intAlongX) cout << "MELALinearInterpFunc(" << GetName() << ")::interpolateFcn integrates using code = " << code << ", coderem = " << coderem << endl;
157  else cout << "MELALinearInterpFunc(" << GetName() << ")::interpolateFcn evaluates using code = " << code << ", coderem = " << coderem << endl;
158  }
159 
160  // Get bins
161  Int_t xbin=-1, xbinmin=-1, xbinmax=-1;
162  MELALinearInterpFunc::T tx=0, txmin=0, txmax=0;
163  if (!intAlongX){ // Case to just compute the value at x
164  if (!testRangeValidity(theXVar)) return 0;
165  xbin = getWhichBin(theXVar);
166  tx = getTVar(theXVar);
167  }
168  else{ // Case to integrate along x
169  MELALinearInterpFunc::T coordmin = theXVar.min(rangeName); cropValueForRange(coordmin);
170  MELALinearInterpFunc::T coordmax = theXVar.max(rangeName); cropValueForRange(coordmax);
171  xbinmin = getWhichBin(coordmin);
172  txmin = getTVar(coordmin);
173  xbinmax = getWhichBin(coordmax);
174  txmax = getTVar(coordmax);
175  }
176 
177  int nxbins = npoints()-1;
178  if (nxbins==0){
180  if (coderem==0) fcnval = dynamic_cast<const RooAbsReal*>(FcnList.at(0))->getVal();
181  else fcnval = dynamic_cast<const RooAbsReal*>(FcnList.at(0))->analyticalIntegral(coderem, rangeName);
182  res = fcnval;
183  if (intAlongX) res *= (txmax-txmin);
184  }
185  else{
186  for (int ix=0; ix<nxbins; ix++){
187  if (
188  (xbin>=0 && ix!=xbin)
189  ||
190  (xbinmin>=0 && xbinmax>=xbinmin && !(xbinmin<=ix && ix<=xbinmax))
191  ) continue;
192 
193  MELALinearInterpFunc::T txlow=0, txhigh=1;
194  if (intAlongX){
195  if (ix==xbinmin) txlow=txmin;
196  if (ix==xbinmax) txhigh=txmax;
197  }
198  else txhigh=tx;
199 
200  MELALinearInterpFunc::T fcnval[2]={ 0 };
201  if (coderem==0){
202  for (unsigned int j=0; j<2; j++) fcnval[j] = dynamic_cast<const RooAbsReal*>(FcnList.at(ix+j))->getVal();
203  }
204  else{
205  for (unsigned int j=0; j<2; j++) fcnval[j] = dynamic_cast<const RooAbsReal*>(FcnList.at(ix+j))->analyticalIntegral(coderem, rangeName);
206  }
207 
209  if (intAlongX) segval += (fcnval[1] - fcnval[0])*(pow(txhigh, 2)-pow(txlow, 2))/2.;
210  else segval += (fcnval[1] - fcnval[0])*txhigh;
211  res += segval;
212 
214  << "MELALinearInterpFunc(" << GetName() << ")::interpolateFcn evaluated bin " << ix
215  << " with txlow = " << txlow << ", txhigh = " << txhigh << ", fcnval[0] = " << fcnval[0] << ", fcnval[1] = " << fcnval[1]
216  << endl;
217  }
218  }
219 
220  return res;
221 }

◆ npoints()

unsigned int MELALinearInterpFunc::npoints ( ) const
inlineprotected

Definition at line 67 of file MELALinearInterpFunc.h.

67 { return XList.size(); }

◆ setEvalFloor()

void MELALinearInterpFunc::setEvalFloor ( MELALinearInterpFunc::T  val)

Definition at line 105 of file MELALinearInterpFunc.cc.

105 { floorEval = val; }

◆ setIntFloor()

void MELALinearInterpFunc::setIntFloor ( MELALinearInterpFunc::T  val)

Definition at line 106 of file MELALinearInterpFunc.cc.

106 { floorInt = val; }

◆ setRangeValidity()

void MELALinearInterpFunc::setRangeValidity ( const T  valmin,
const T  valmax 
)

Definition at line 306 of file MELALinearInterpFunc.cc.

306  {
307  T* range[2];
308  range[0] = &rangeXmin;
309  range[1] = &rangeXmax;
310  *(range[0])=valmin;
311  *(range[1])=valmax;
312 }

◆ setVerbosity()

void MELALinearInterpFunc::setVerbosity ( VerbosityLevel  flag)

Definition at line 104 of file MELALinearInterpFunc.cc.

104 { verbosity = flag; }

◆ testRangeValidity()

Bool_t MELALinearInterpFunc::testRangeValidity ( const T val) const
protected

Definition at line 300 of file MELALinearInterpFunc.cc.

300  {
301  const T* range[2];
302  range[0] = &rangeXmin;
303  range[1] = &rangeXmax;
304  return (*(range[0])>*(range[1]) || (val>=*(range[0]) && val<=*(range[1])));
305 }

Member Data Documentation

◆ FcnList

RooListProxy MELALinearInterpFunc::FcnList
protected

Definition at line 64 of file MELALinearInterpFunc.h.

◆ floorEval

T MELALinearInterpFunc::floorEval
protected

Definition at line 56 of file MELALinearInterpFunc.h.

◆ floorInt

T MELALinearInterpFunc::floorInt
protected

Definition at line 57 of file MELALinearInterpFunc.h.

◆ leafDepsList

RooListProxy MELALinearInterpFunc::leafDepsList
protected

Definition at line 65 of file MELALinearInterpFunc.h.

◆ rangeXmax

T MELALinearInterpFunc::rangeXmax
protected

Definition at line 60 of file MELALinearInterpFunc.h.

◆ rangeXmin

T MELALinearInterpFunc::rangeXmin
protected

Definition at line 59 of file MELALinearInterpFunc.h.

◆ theXVar

RooRealProxy MELALinearInterpFunc::theXVar
protected

Definition at line 63 of file MELALinearInterpFunc.h.

◆ useFloor

Bool_t MELALinearInterpFunc::useFloor
protected

Definition at line 55 of file MELALinearInterpFunc.h.

◆ verbosity

VerbosityLevel MELALinearInterpFunc::verbosity
protected

Definition at line 54 of file MELALinearInterpFunc.h.

◆ XList

std::vector<T> MELALinearInterpFunc::XList
protected

Definition at line 62 of file MELALinearInterpFunc.h.


The documentation for this class was generated from the following files:
MELALinearInterpFunc::kVerbose
@ kVerbose
Definition: MELALinearInterpFunc.h:20
MELALinearInterpFunc::FcnList
RooListProxy FcnList
Definition: MELALinearInterpFunc.h:64
MELALinearInterpFunc::verbosity
VerbosityLevel verbosity
Definition: MELALinearInterpFunc.h:54
value
pymela::gHIGGS_KAPPA value("gHIGGS_KAPPA_TILDE", pymela::gHIGGS_KAPPA_TILDE) .value("SIZE_HQQ"
TNumericUtil::KahanAccumulator
Definition: MELAAccumulators.h:27
MELALinearInterpFunc::leafDepsList
RooListProxy leafDepsList
Definition: MELALinearInterpFunc.h:65
MELALinearInterpFunc::getWhichBin
Int_t getWhichBin(const T &val) const
Definition: MELALinearInterpFunc.cc:109
MELALinearInterpFunc::cropValueForRange
void cropValueForRange(T &val) const
Definition: MELALinearInterpFunc.cc:313
MELALinearInterpFunc::rangeXmax
T rangeXmax
Definition: MELALinearInterpFunc.h:60
MELALinearInterpFunc::MELALinearInterpFunc
MELALinearInterpFunc()
MELALinearInterpFunc::kSilent
@ kSilent
Definition: MELALinearInterpFunc.h:18
modparameters::kappa
complex(8), public kappa
Definition: mod_Parameters.F90:882
MELALinearInterpFunc::getTVar
T getTVar(const T &val) const
Definition: MELALinearInterpFunc.cc:140
MELALinearInterpFunc::theXVar
RooRealProxy theXVar
Definition: MELALinearInterpFunc.h:63
MELALinearInterpFunc::npoints
unsigned int npoints() const
Definition: MELALinearInterpFunc.h:67
MELALinearInterpFunc::floorEval
T floorEval
Definition: MELALinearInterpFunc.h:56
MELALinearInterpFunc::rangeXmin
T rangeXmin
Definition: MELALinearInterpFunc.h:59
MELALinearInterpFunc::T
Float_t T
Definition: MELALinearInterpFunc.h:15
dd_global::cout
integer cout
Definition: DD_global.F90:21
MELALinearInterpFunc::floorInt
T floorInt
Definition: MELALinearInterpFunc.h:57
MELALinearInterpFunc::XList
std::vector< T > XList
Definition: MELALinearInterpFunc.h:62
MELALinearInterpFunc::analyticalIntegral
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Definition: MELALinearInterpFunc.cc:290
MELALinearInterpFunc::useFloor
Bool_t useFloor
Definition: MELALinearInterpFunc.h:55
MELALinearInterpFunc::kError
@ kError
Definition: MELALinearInterpFunc.h:19
MELALinearInterpFunc::getKappa
T getKappa(const Int_t &bin) const
Definition: MELALinearInterpFunc.cc:130
MELALinearInterpFunc::interpolateFcn
T interpolateFcn(Int_t code, const char *rangeName=0) const
Definition: MELALinearInterpFunc.cc:147
MELALinearInterpFunc::testRangeValidity
Bool_t testRangeValidity(const T &val) const
Definition: MELALinearInterpFunc.cc:300