REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestComponent.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see https://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see https://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
41#include "TRestComponent.h"
42
43#include <TKey.h>
44#include <TLatex.h>
45
46#include <numeric>
47
48ClassImp(TRestComponent);
49
54
69TRestComponent::TRestComponent(const char* cfgFileName, const std::string& name)
70 : TRestMetadata(cfgFileName) {
71 RESTDebug << "Entering TRestComponent constructor( cfgFileName, name )" << RESTendl;
72 RESTDebug << "File: " << cfgFileName << " Name: " << name << RESTendl;
73}
74
79
87 // SetSectionName(this->ClassName());
88
90 if (!fNodeDensity.empty() && fRandom) return;
91
92 if (!fRandom) {
93 delete fRandom;
94 fRandom = nullptr;
95 }
96
97 fRandom = new TRandom3(fSeed);
98 fSeed = fRandom->TRandom::GetSeed();
99
100 if (fStepParameterValue > 0) {
103 } else {
104 if (!fParameterizationNodes.empty()) FillHistograms();
105 }
106}
107
115 fNodeDensity.clear();
116
117 fSeed = seed;
118 FillHistograms();
119}
120
126void TRestComponent::RegenerateParametricNodes(Double_t from, Double_t to, Double_t step,
127 Bool_t expIncrease) {
128 fStepParameterValue = step;
131 fExponential = expIncrease;
132
134
135 if (expIncrease) {
136 for (double p = from; p < to; p *= step) fParameterizationNodes.push_back(p);
137 } else {
138 for (double p = from; p < to; p += step) fParameterizationNodes.push_back(p);
139 }
140
141 if (fParameterizationNodes.empty()) return;
143}
144
149Int_t TRestComponent::GetVariableIndex(std::string varName) {
150 int n = 0;
151 for (const auto& v : fVariables) {
152 if (v == varName) return n;
153 n++;
154 }
155
156 return -1;
157}
158
168Double_t TRestComponent::GetRate(std::vector<Double_t> point) {
169 if (!fResponse) {
170 return GetRawRate(point);
171 }
172
173 std::string responseVariable = fResponse->GetVariable();
174 Int_t respVarIndex = GetVariableIndex(responseVariable);
175
176 if (respVarIndex == -1) {
177 RESTError << "The response variable `" << responseVariable << "`, defined inside TRestResponse,"
178 << RESTendl;
179 RESTError << "could not be found inside the component." << RESTendl;
180 RESTError << "Please, check the component variable names." << RESTendl;
181 return 0;
182 }
183
184 std::vector<std::pair<Double_t, Double_t> > response = fResponse->GetResponse(point[respVarIndex]);
185
186 Double_t rate = 0;
187 for (const auto& resp : response) {
188 std::vector<Double_t> newPoint = point;
189 newPoint[respVarIndex] = resp.first;
190 rate += resp.second * GetRawRate(newPoint);
191 }
192
193 return rate;
194}
195
217Double_t TRestComponent::GetNormalizedRate(std::vector<Double_t> point) {
218 Double_t normFactor = 1;
219 for (size_t n = 0; n < GetDimensions(); n++) normFactor *= fNbins[n] / (fRanges[n].Y() - fRanges[n].X());
220
221 return normFactor * GetRate(point);
222}
223
242Double_t TRestComponent::GetRawRate(std::vector<Double_t> point) {
243 if (point.size() != GetDimensions()) {
244 RESTError << "The size of the point given is : " << point.size() << RESTendl;
245 RESTError << "The density distribution dimensions are : " << GetDimensions() << RESTendl;
246 RESTError << "Point must have the same dimension as the distribution" << RESTendl;
247 return 0;
248 }
249
250 if (!HasNodes() && !HasDensity()) {
251 RESTError << "TRestComponent::GetRawRate. The component has no nodes!" << RESTendl;
252 RESTError << "Try calling TRestComponent::Initialize" << RESTendl;
253
254 RESTInfo << "Trying to initialize" << RESTendl;
255 Initialize();
256 if (HasNodes())
257 RESTInfo << "Sucess!" << RESTendl;
258 else
259 return 0;
260 }
261
262 if (HasNodes() && fActiveNode == -1) {
263 RESTError << "TRestComponent::GetRawRate. Active node has not been defined" << RESTendl;
264 return 0;
265 }
266
267 for (size_t n = 0; n < point.size(); n++) {
268 // The point is outside boundaries
269 if (point[n] < fRanges[n].X() || point[n] > fRanges[n].Y()) return 0;
270 }
271
272 Int_t centerBin[GetDimensions()];
273 Double_t centralDensity = GetDensity()->GetBinContent(GetDensity()->GetBin(point.data()), centerBin);
274 if (!Interpolation()) return centralDensity;
275
276 std::vector<Int_t> direction;
277 std::vector<Double_t> nDist;
278 for (size_t dim = 0; dim < GetDimensions(); dim++) {
279 Double_t x1 = GetBinCenter(dim, centerBin[dim] - 1);
280 Double_t x2 = GetBinCenter(dim, centerBin[dim] + 1);
281
282 if (centerBin[dim] == 1 || centerBin[dim] == fNbins[dim]) {
283 direction.push_back(0);
284 nDist.push_back(0);
285 } else if (x2 - point[dim] > point[dim] - x1) {
286 // we chose left bin (x1) since it is closer than right bin
287 direction.push_back(-1);
288 nDist.push_back(1 - 2 * (point[dim] - x1) / (x2 - x1));
289 } else {
290 direction.push_back(1);
291 nDist.push_back(1 - 2 * (x2 - point[dim]) / (x2 - x1));
292 }
293 }
294
295 // In 3-dimensions we got 8 points to interpolate
296 // In 4-dimensions we would get 16 points to interpolate
297 // ...
298 Int_t nPoints = (Int_t)TMath::Power(2, (Int_t)GetDimensions());
299
300 Double_t sum = 0;
301 for (int n = 0; n < nPoints; n++) {
302 std::vector<int> cell = REST_StringHelper::IntegerToBinary(n, GetDimensions());
303
304 Double_t weightDistance = 1;
305 int cont = 0;
306 for (const auto& c : cell) {
307 if (c == 0)
308 weightDistance *= (1 - nDist[cont]);
309 else
310 weightDistance *= nDist[cont];
311 cont++;
312 }
313
314 for (size_t k = 0; k < cell.size(); k++) cell[k] = cell[k] * direction[k] + centerBin[k];
315
316 Double_t density = GetDensity()->GetBinContent(cell.data());
317 sum += density * weightDistance;
318 }
319
320 return sum;
321}
322
328 THnD* dHist = GetDensityForActiveNode();
329 if (!dHist) return 0;
330
331 Double_t integral = 0;
332 for (Int_t n = 0; n < dHist->GetNbins(); ++n) {
333 Int_t centerBin[GetDimensions()];
334 std::vector<Double_t> point;
335
336 dHist->GetBinContent(n, centerBin);
337 for (size_t d = 0; d < GetDimensions(); ++d) point.push_back(GetBinCenter(d, centerBin[d]));
338
339 Bool_t skip = false;
340 for (size_t d = 0; d < GetDimensions(); ++d) {
341 if (point[d] < fRanges[d].X() || point[d] > fRanges[d].Y()) skip = true;
342 }
343 if (!skip) integral += GetRate(point);
344 }
345
346 return integral;
347}
348
354 Double_t maxRate = 0;
355 for (size_t n = 0; n < fParameterizationNodes.size(); n++) {
356 SetActiveNode((Int_t)n);
357 Double_t rate = GetTotalRate();
358 if (rate > maxRate) maxRate = rate;
359 }
360 return maxRate;
361}
362
368 Double_t rate = 0;
369 for (size_t n = 0; n < fParameterizationNodes.size(); n++) {
370 SetActiveNode((Int_t)n);
371 rate += GetTotalRate();
372 }
373 return rate;
374}
375
381Double_t TRestComponent::GetBinCenter(Int_t nDim, const Int_t bin) {
382 return fRanges[nDim].X() + (fRanges[nDim].Y() - fRanges[nDim].X()) * ((double)bin - 0.5) / fNbins[nDim];
383}
384
385ROOT::RVecD TRestComponent::GetRandom() {
386 Double_t* tuple = new Double_t[GetDimensions()];
387
388 ROOT::RVecD result;
389 if (!GetDensity()) {
390 for (size_t n = 0; n < GetDimensions(); n++) result.push_back(0);
391 RESTWarning << "TRestComponent::GetRandom. Component might not be initialized! Use "
392 "TRestComponent::Initialize()."
393 << RESTendl;
394 return result;
395 }
396
397 GetDensity()->GetRandom(tuple);
398
399 for (size_t n = 0; n < GetDimensions(); n++) result.push_back(tuple[n]);
400 return result;
401}
402
403ROOT::RDF::RNode TRestComponent::GetMonteCarloDataFrame(Int_t N) {
404 ROOT::RDF::RNode df = ROOT::RDataFrame(N);
405
406 // Function to fill the RDataFrame using GetRandom method
407 auto fillRndm = [&]() {
408 ROOT::RVecD randomValues = GetRandom();
409 return randomValues;
410 };
411 df = df.Define("Rndm", fillRndm);
412
413 // Creating dedicated columns for each GetRandom component
414 for (size_t i = 0; i < fVariables.size(); ++i) {
415 auto varName = fVariables[i];
416 auto FillRand = [i](const ROOT::RVecD& randomValues) { return randomValues[i]; };
417 df = df.Define(varName, FillRand, {"Rndm"});
418 }
419
420 /* Excluding Rndm from df */
421 std::vector<std::string> columns = df.GetColumnNames();
422 std::vector<std::string> excludeColumns = {"Rndm"};
423 columns.erase(std::remove_if(columns.begin(), columns.end(),
424 [&excludeColumns](std::string elem) {
425 return std::find(excludeColumns.begin(), excludeColumns.end(), elem) !=
426 excludeColumns.end();
427 }),
428 columns.end());
429
430 std::string user = getenv("USER");
431 std::string fOutName = "/tmp/rest_output_" + user + ".root";
432 df.Snapshot("AnalysisTree", fOutName, columns);
433
434 df = ROOT::RDataFrame("AnalysisTree", fOutName);
435
436 return df;
437}
438
448TCanvas* TRestComponent::DrawComponent(std::vector<std::string> drawVariables,
449 std::vector<std::string> scanVariables, Int_t binScanSize,
450 TString drawOption) {
451 if (drawVariables.size() > 2 || drawVariables.size() == 0) {
452 RESTError << "TRestComponent::DrawComponent. The number of variables to be drawn must "
453 "be 1 or 2!"
454 << RESTendl;
455 return fCanvas;
456 }
457
458 if (scanVariables.size() > 2 || scanVariables.size() == 0) {
459 RESTError << "TRestComponent::DrawComponent. The number of variables to be scanned must "
460 "be 1 or 2!"
461 << RESTendl;
462 return fCanvas;
463 }
464
466 std::vector<int> scanIndexes;
467 for (const auto& x : scanVariables) scanIndexes.push_back(GetVariableIndex(x));
468
469 Int_t nPlots = 1;
470 size_t n = 0;
471 for (const auto& x : scanIndexes) {
472 if (fNbins[x] % binScanSize != 0) {
473 RESTWarning << "The variable " << scanVariables[n] << " contains " << fNbins[x]
474 << " bins and it doesnt match with a bin size " << binScanSize << RESTendl;
475 RESTWarning << "The bin size must be a multiple of the number of bins." << RESTendl;
476 RESTWarning << "Redefining bin size to 1." << RESTendl;
477 binScanSize = 1;
478 }
479 nPlots *= fNbins[x] / binScanSize;
480 n++;
481 }
482
484 Int_t nPlotsX = 0;
485 Int_t nPlotsY = 0;
486
487 if (scanIndexes.size() == 2) {
488 nPlotsX = fNbins[scanIndexes[0]] / binScanSize;
489 nPlotsY = fNbins[scanIndexes[1]] / binScanSize;
490 } else {
491 nPlotsX = TRestTools::CanvasDivisions(nPlots)[1];
492 nPlotsY = TRestTools::CanvasDivisions(nPlots)[0];
493 }
494
495 RESTInfo << "Number of plots to be generated: " << nPlots << RESTendl;
496 RESTInfo << "Canvas size : " << nPlotsX << " x " << nPlotsY << RESTendl;
497
499 if (fCanvas != nullptr) {
500 delete fCanvas;
501 fCanvas = nullptr;
502 }
503
504 fCanvas = new TCanvas(this->GetName(), this->GetName(), 0, 0, nPlotsX * 640, nPlotsY * 480);
505 fCanvas->Divide(nPlotsX, nPlotsY, 0.01, 0.01);
506
507 std::vector<int> variableIndexes;
508 for (const auto& x : drawVariables) variableIndexes.push_back(GetVariableIndex(x));
509
510 for (int n = 0; n < nPlotsX; n++)
511 for (int m = 0; m < nPlotsY; m++) {
512 TPad* pad = (TPad*)fCanvas->cd(n * nPlotsY + m + 1);
513 pad->SetFixedAspectRatio(true);
514
515 THnD* hnd = GetDensity();
516
517 int binXo = binScanSize * n + 1;
518 int binXf = binScanSize * n + binScanSize;
519 int binYo = binScanSize * m + 1;
520 int binYf = binScanSize * m + binScanSize;
521
522 if (scanVariables.size() == 2) {
523 hnd->GetAxis(scanIndexes[0])->SetRange(binXo, binXf);
524 hnd->GetAxis(scanIndexes[1])->SetRange(binYo, binYf);
525 } else if (scanVariables.size() == 1) {
526 binXo = binScanSize * nPlotsY * n + binScanSize * m + 1;
527 binXf = binScanSize * nPlotsY * n + binScanSize * m + binScanSize;
528 hnd->GetAxis(scanIndexes[0])->SetRange(binXo, binXf);
529 }
530
531 if (variableIndexes.size() == 1) {
532 TH1D* h1 = hnd->Projection(variableIndexes[0]);
533 std::string hName;
534
535 if (scanIndexes.size() == 2)
536 hName = scanVariables[0] + "(" + IntegerToString(binXo) + ", " + IntegerToString(binXf) +
537 ") " + scanVariables[1] + "(" + IntegerToString(binYo) + ", " +
538 IntegerToString(binYf) + ") ";
539
540 if (scanIndexes.size() == 1)
541 hName = scanVariables[0] + "(" + IntegerToString(binXo) + ", " + IntegerToString(binXf) +
542 ") ";
543
544 TH1D* newh = (TH1D*)h1->Clone(hName.c_str());
545 newh->SetTitle(hName.c_str());
546 newh->SetStats(false);
547 newh->GetXaxis()->SetTitle((TString)drawVariables[0]);
548 newh->SetMarkerStyle(kFullCircle);
549 newh->Draw("PLC PMC");
550
551 TString entriesStr = "Entries: " + IntegerToString(newh->GetEntries());
552 TLatex* textLatex = new TLatex(0.62, 0.825, entriesStr);
553 textLatex->SetNDC();
554 textLatex->SetTextColor(1);
555 textLatex->SetTextSize(0.05);
556 textLatex->Draw("same");
557 delete h1;
558 }
559
560 if (variableIndexes.size() == 2) {
561 TH2D* h2 = hnd->Projection(variableIndexes[0], variableIndexes[1]);
562
563 std::string hName;
564 if (scanIndexes.size() == 2)
565 hName = scanVariables[0] + "(" + IntegerToString(binXo) + ", " + IntegerToString(binXf) +
566 ") " + scanVariables[1] + "(" + IntegerToString(binYo) + ", " +
567 IntegerToString(binYf) + ") ";
568
569 if (scanIndexes.size() == 1)
570 hName = scanVariables[0] + "(" + IntegerToString(binXo) + ", " + IntegerToString(binXf) +
571 ") ";
572
573 TH2D* newh = (TH2D*)h2->Clone(hName.c_str());
574 newh->SetStats(false);
575 newh->GetXaxis()->SetTitle((TString)drawVariables[0]);
576 newh->GetYaxis()->SetTitle((TString)drawVariables[1]);
577 newh->SetTitle(hName.c_str());
578 newh->Draw(drawOption);
579
580 TString entriesStr = "Entries: " + IntegerToString(newh->GetEntries());
581 TLatex* textLatex = new TLatex(0.62, 0.825, entriesStr);
582 textLatex->SetNDC();
583 textLatex->SetTextColor(1);
584 textLatex->SetTextSize(0.05);
585 textLatex->Draw("same");
586 delete h2;
587 }
588 }
589
590 return fCanvas;
591}
592
596void TRestComponent::LoadResponse(const TRestResponse& resp) {
597 if (fResponse) {
598 delete fResponse;
599 fResponse = nullptr;
600 }
601
602 fResponse = (TRestResponse*)resp.Clone("response");
604
606}
607
613
614 RESTMetadata << "Component nature : " << fNature << RESTendl;
615 RESTMetadata << " " << RESTendl;
616
617 RESTMetadata << "Random seed : " << fSeed << RESTendl;
618 if (fSamples) RESTMetadata << "Samples : " << fSamples << RESTendl;
619 RESTMetadata << " " << RESTendl;
620
621 if (fVariables.size() != fRanges.size())
622 RESTWarning << "The number of variables does not match with the number of defined ranges!"
623 << RESTendl;
624
625 else if (fVariables.size() != fNbins.size())
626 RESTWarning << "The number of variables does not match with the number of defined bins!" << RESTendl;
627 else {
628 int n = 0;
629 RESTMetadata << " === Variables === " << RESTendl;
630 for (const auto& varName : fVariables) {
631 RESTMetadata << " - Name: " << varName << " Range: (" << fRanges[n].X() << ", " << fRanges[n].Y()
632 << ") bins: " << fNbins[n] << RESTendl;
633 n++;
634 }
635 }
636
637 if (!fParameterizationNodes.empty()) {
638 RESTMetadata << " " << RESTendl;
639 RESTMetadata << " === Parameterization === " << RESTendl;
640 RESTMetadata << "- Parameter : " << fParameter << RESTendl;
641
642 RESTMetadata << " - Number of parametric nodes : " << fParameterizationNodes.size() << RESTendl;
643 RESTMetadata << " " << RESTendl;
644 RESTMetadata << " Use : PrintNodes() for additional info" << RESTendl;
645
646 if (fStepParameterValue > 0) {
647 RESTMetadata << " " << RESTendl;
648 RESTMetadata << " Nodes were automatically generated using these parameters" << RESTendl;
649 RESTMetadata << " - First node : " << fFirstParameterValue << RESTendl;
650 RESTMetadata << " - Upper limit node : " << fLastParameterValue << RESTendl;
651 RESTMetadata << " - Increasing step : " << fStepParameterValue << RESTendl;
652 if (fExponential)
653 RESTMetadata << " - Increases exponentially" << RESTendl;
654 else
655 RESTMetadata << " - Increases linearly" << RESTendl;
656 }
657 }
658
659 if (fResponse) {
660 RESTMetadata << " " << RESTendl;
661 RESTMetadata << "A response matrix was loaded inside the component" << RESTendl;
662 RESTMetadata << "You may get more details using TRestComponent::GetResponse()->PrintMetadata()"
663 << RESTendl;
664 }
665
666 RESTMetadata << "----" << RESTendl;
667}
668
673 std::cout << " - Number of nodes : " << fParameterizationNodes.size() << std::endl;
674 for (const auto& x : fParameterizationNodes) std::cout << x << " ";
675 std::cout << std::endl;
676}
677
683
684 auto ele = GetElement("cVariable");
685 while (ele != nullptr) {
686 std::string name = GetParameter("name", ele, "");
687 TVector2 v = Get2DVectorParameterWithUnits("range", ele, TVector2(-1, -1));
688 Int_t bins = StringToInteger(GetParameter("bins", ele, "0"));
689
690 if (name.empty() || (v.X() == -1 && v.Y() == -1) || bins == 0) {
691 RESTWarning << "TRestComponentFormula::fVariable. Problem with definition." << RESTendl;
692 RESTWarning << "Name: " << name << " range: (" << v.X() << ", " << v.Y() << ") bins: " << bins
693 << RESTendl;
694 } else {
695 fVariables.push_back(name);
696 fRanges.push_back(v);
697 fNbins.push_back(bins);
698 }
699
700 ele = GetNextElement(ele);
701 }
702
703 if (fNbins.size() == 0)
704 RESTError << "TRestComponent::InitFromConfigFile. No cVariables where found!" << RESTendl;
705
706 if (fResponse) {
707 delete fResponse;
708 fResponse = nullptr;
709 }
710
713}
714
719Int_t TRestComponent::FindActiveNode(Double_t node) {
720 int n = 0;
721 for (const auto& v : fParameterizationNodes) {
722 Double_t pUp = node * (1 + fPrecision / 2);
723 Double_t pDown = node * (1 - fPrecision / 2);
724 if (v > pDown && v < pUp) {
725 fActiveNode = n;
726 return fActiveNode;
727 }
728 n++;
729 }
730
731 return -1;
732}
733
738Int_t TRestComponent::SetActiveNode(Double_t node) {
739 int n = 0;
740 for (const auto& v : fParameterizationNodes) {
741 Double_t pUp = node * (1 + fPrecision / 2);
742 Double_t pDown = node * (1 - fPrecision / 2);
743 if (v > pDown && v < pUp) {
744 fActiveNode = n;
745 return fActiveNode;
746 }
747 n++;
748 }
749
750 RESTError << "Parametric node : " << node << " was not found in component" << RESTendl;
751 RESTError << "Keeping previous node as active : " << fParameterizationNodes[fActiveNode] << RESTendl;
752 PrintNodes();
753
754 return fActiveNode;
755}
756
760THnD* TRestComponent::GetDensityForNode(Double_t node) {
761 int n = 0;
762 for (const auto& x : fParameterizationNodes) {
763 if (x == node) {
764 return fNodeDensity[n];
765 }
766 n++;
767 }
768
769 RESTError << "Parametric node : " << node << " was not found in component" << RESTendl;
770 PrintNodes();
771 return nullptr;
772}
773
777THnD* TRestComponent::GetDensityForActiveNode() {
778 if (fActiveNode >= 0) return fNodeDensity[fActiveNode];
779
780 RESTError << "The active node is invalid" << RESTendl;
781 PrintNodes();
782 return nullptr;
783}
784
789TH1D* TRestComponent::GetHistogram(Double_t node, std::string varName) {
790 SetActiveNode(node);
791 return GetHistogram(varName);
792}
793
799TH1D* TRestComponent::GetHistogram(std::string varName) {
800 if (fActiveNode < 0) return nullptr;
801
802 Int_t v1 = GetVariableIndex(varName);
803
804 if (v1 >= 0 && GetDensityForActiveNode()) return GetDensityForActiveNode()->Projection(v1);
805
806 return nullptr;
807}
808
813TH2D* TRestComponent::GetHistogram(Double_t node, std::string varName1, std::string varName2) {
814 SetActiveNode(node);
815 return GetHistogram(varName1, varName2);
816}
817
823TH2D* TRestComponent::GetHistogram(std::string varName1, std::string varName2) {
824 if (fActiveNode < 0) return nullptr;
825
826 Int_t v1 = GetVariableIndex(varName1);
827 Int_t v2 = GetVariableIndex(varName2);
828
829 if (v1 >= 0 && v2 >= 0)
830 if (GetDensityForActiveNode()) return GetDensityForActiveNode()->Projection(v1, v2);
831
832 return nullptr;
833}
834
839TH3D* TRestComponent::GetHistogram(Double_t node, std::string varName1, std::string varName2,
840 std::string varName3) {
841 SetActiveNode(node);
842 return GetHistogram(varName1, varName2, varName3);
843}
844
850TH3D* TRestComponent::GetHistogram(std::string varName1, std::string varName2, std::string varName3) {
851 if (fActiveNode < 0) return nullptr;
852
853 Int_t v1 = GetVariableIndex(varName1);
854 Int_t v2 = GetVariableIndex(varName2);
855 Int_t v3 = GetVariableIndex(varName3);
856
857 if (v1 >= 0 && v2 >= 0 && v3 >= 0)
858 if (GetDensityForActiveNode()) return GetDensityForActiveNode()->Projection(v1, v2, v3);
859
860 return nullptr;
861}
It defines a background/signal model distribution in a given parameter space (tipically x,...
Double_t fFirstParameterValue
It defines the first parametric node value in case of automatic parameter generation.
Double_t fLastParameterValue
It defines the upper limit for the automatic parametric node values generation.
UInt_t fSeed
Seed used in random generator.
Double_t fStepParameterValue
It defines the increasing step for automatic parameter list generation.
void InitFromConfigFile() override
It customizes the retrieval of XML data values of this class.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
Double_t GetAllNodesIntegratedRate()
This method returns the integrated total rate for all the nodes The result will be returned in s-1.
Int_t fActiveNode
It is used to define the node that will be accessed for rate retrieval.
Int_t FindActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
void RegenerateParametricNodes(Double_t from, Double_t to, Double_t step, Bool_t expIncrease=false)
It allows to produce a parameter nodes list providing the initial value, the final value and the step...
Int_t fSamples
It introduces a fixed number of samples (if 0 it will take all available samples)
Int_t GetVariableIndex(std::string varName)
It returns the position of the fVariable element for the variable name given by argument.
Double_t GetBinCenter(Int_t nDim, const Int_t bin)
It returns the bin center of the given component dimension.
void PrintNodes()
It prints out on screen the values of the parametric node.
Double_t GetRawRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
Bool_t fExponential
It true the parametric values automatically generated will grow exponentially.
std::string fParameter
It is used to parameterize a set of distribution densities (e.g. WIMP or axion mass)
Double_t GetNormalizedRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
TRestComponent()
Default constructor.
TH1D * GetHistogram(Double_t node, std::string varName)
It returns a 1-dimensional projected histogram for the variable names provided in the argument.
Float_t fPrecision
A precision used to select the node value with a given range defined as a fraction of the value.
std::vector< Int_t > fNbins
The number of bins in which we should divide each variable.
Double_t GetMaxRate()
This method returns the total rate for the node that has the highest contribution The result will be ...
TCanvas * DrawComponent(std::vector< std::string > drawVariables, std::vector< std::string > scanVariables, Int_t binScanSize=1, TString drawOption="")
A method allowing to draw a series of plots representing the density distributions.
void Initialize() override
It initializes the random number. We avoid to define the section name here since we will never define...
TRestResponse * fResponse
A pointer to the detector response.
Double_t GetTotalRate()
This method integrates the rate to all the parameter space defined in the density function....
std::vector< TVector2 > fRanges
The range of each of the variables used to create the PDF distribution.
std::vector< Double_t > fParameterizationNodes
It defines the nodes of the parameterization (Initialized by the dataset)
Bool_t HasNodes()
It returns true if any nodes have been defined.
TCanvas * fCanvas
A canvas for drawing the active node component.
Double_t GetRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
TRandom3 * fRandom
Internal process random generator.
std::vector< std::string > fVariables
A list with the branches that will be used to create the distribution space.
~TRestComponent()
Default destructor.
void RegenerateHistograms(UInt_t seed=0)
It will produce a histogram with the distribution defined using the variables and the weights for eac...
Int_t SetActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
std::vector< THnD * > fNodeDensity
The generated N-dimensional variable space density for a given node.
std::string fNature
It defines the component type (unknown/signal/background)
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
TiXmlElement * GetElement(std::string eleDeclare, TiXmlElement *e=nullptr)
Get an xml element from a given parent element, according to its declaration.
TRestMetadata * InstantiateChildMetadata(int index, std::string pattern="")
This method will retrieve a new TRestMetadata instance of a child element of the present TRestMetadat...
virtual void InitFromConfigFile()
To make settings from rml file. This method must be implemented in the derived class.
TiXmlElement * GetNextElement(TiXmlElement *e)
Get the next sibling xml element of this element, with same eleDeclare.
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
A response matrix that might be applied to a given component inside a TRestComponent.
Definition: TRestResponse.h:29
void LoadResponse(Bool_t transpose=true)
It loads into the fResponseMatrix data member the response from a file.
std::vector< std::pair< Double_t, Double_t > > GetResponse(Double_t input)
This method will return a vector of std::pair, each pair will contain the output energy together with...
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
static std::vector< int > CanvasDivisions(int n)
It returns a vector with 2 components {a,b}, the components satisfy that a x b = n,...
std::string IntegerToString(Int_t n, std::string format="%d")
Gets a string from an integer.
Int_t StringToInteger(std::string in)
Gets an integer from a string.
std::vector< int > IntegerToBinary(int number, size_t dimension=0)
It returns an integer vector with the binary digits decomposed.