REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionSolarQCDFlux.cxx
1/******************** REST disclaimer ***********************************
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 http://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 http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
211
212#include "TRestAxionSolarQCDFlux.h"
213using namespace std;
214
215ClassImp(TRestAxionSolarQCDFlux);
216
221
236TRestAxionSolarQCDFlux::TRestAxionSolarQCDFlux(const char* cfgFileName, string name)
237 : TRestAxionSolarFlux(cfgFileName) {
239
241}
242
247
253 if (fFluxDataFile == "" && fFluxSptFile == "") return false;
254
255 if (fTablesLoaded) return true;
256
258 ReadFluxFile();
259 } else {
262 }
263
265
266 fTablesLoaded = true;
267
268 return true;
269}
270
276 if (fFluxDataFile == "") {
277 RESTDebug << "TRestAxionSolarflux::LoadContinuumFluxTable. No solar flux table was defined"
278 << RESTendl;
279 return;
280 }
281
282 string fullPathName = SearchFile((string)fFluxDataFile);
283
284 RESTDebug << "Loading table from file : " << RESTendl;
285 RESTDebug << "File : " << fullPathName << RESTendl;
286
287 std::vector<std::vector<Float_t>> fluxTable;
289 std::vector<std::vector<Double_t>> doubleTable;
290 if (!TRestTools::ReadASCIITable(fullPathName, doubleTable)) {
291 RESTError << "TRestAxionSolarQCDFlux::LoadContinuumFluxTable. " << RESTendl;
292 RESTError << "Could not read solar flux table : " << fFluxDataFile << RESTendl;
293 return;
294 }
295 for (const auto& row : doubleTable) {
296 std::vector<Float_t> floatVec(row.begin(), row.end());
297 fluxTable.push_back(floatVec);
298 }
300 TRestTools::ReadBinaryTable(fullPathName, fluxTable);
301 else {
302 fluxTable.clear();
303 RESTError << "Filename extension was not recognized!" << RESTendl;
304 RESTError << "Solar flux table will not be populated" << RESTendl;
305 RESTError << "Filename extension: " << TRestTools::GetFileNameExtension(fFluxDataFile) << RESTendl;
306 }
307
308 if (fluxTable.size() != 100 && fluxTable[0].size() != 200) {
309 fluxTable.clear();
310 RESTError << "LoadContinuumFluxTable. The table does not contain the right number of rows or columns"
311 << RESTendl;
312 RESTError << "Table will not be populated" << RESTendl;
313 }
314
315 for (unsigned int n = 0; n < fluxTable.size(); n++) {
316 TH1F* h = new TH1F(Form("%s_ContinuumFluxAtRadius%d", GetName(), n), "", 200, 0, 20);
317 for (unsigned int m = 0; m < fluxTable[n].size(); m++) h->SetBinContent(m + 1, fluxTable[n][m]);
318 fFluxTable.push_back(h);
319 }
320}
321
327 if (fFluxSptFile == "") {
328 RESTDebug << "TRestAxionSolarflux::LoadMonoChromaticFluxTable. No solar flux monochromatic table was "
329 "defined"
330 << RESTendl;
331 return;
332 }
333
334 string fullPathName = SearchFile((string)fFluxSptFile);
335
336 RESTDebug << "Loading monochromatic lines from file : " << RESTendl;
337 RESTDebug << "File : " << fullPathName << RESTendl;
338
339 std::vector<std::vector<Double_t>> doubleTable;
340 if (!TRestTools::ReadASCIITable(fullPathName, doubleTable)) {
341 RESTError << "TRestAxionSolarQCDFlux::LoadMonoChromaticFluxTable." << RESTendl;
342 RESTError << "Could not read solar flux table : " << fFluxSptFile << RESTendl;
343 return;
344 }
345
346 std::vector<std::vector<Float_t>> asciiTable;
347 for (const auto& row : doubleTable) {
348 std::vector<Float_t> floatVec(row.begin(), row.end());
349 asciiTable.push_back(floatVec);
350 }
351
352 fFluxLines.clear();
353
354 if (asciiTable.size() != 101) {
355 RESTError << "LoadMonoChromaticFluxTable. The table does not contain the right number of rows"
356 << RESTendl;
357 RESTError << "Table will not be populated" << RESTendl;
358 return;
359 }
360
361 for (unsigned int en = 0; en < asciiTable[0].size(); en++) {
362 Float_t energy = asciiTable[0][en];
363 TH1F* h = new TH1F(Form("%s_MonochromeFluxAtEnergy%6.4lf", GetName(), energy), "", 100, 0, 1);
364 for (unsigned int r = 1; r < asciiTable.size(); r++) h->SetBinContent(r, asciiTable[r][en]);
365 fFluxLines[energy] = h;
366 }
367}
368
374 if (fBinSize <= 0) {
375 RESTError << "TRestAxionSolarflux::ReadFluxFile. Energy bin size of .flux file must be specified."
376 << RESTendl;
377 RESTError << "Please, define binSize parameter in eV." << RESTendl;
378 return;
379 }
380
381 if (fPeakSigma <= 0) {
382 RESTWarning << "TRestAxionSolarflux::ReadFluxFile. Peak sigma must be specified to generate "
383 "monochromatic spectrum."
384 << RESTendl;
385 RESTWarning << "Only continuum table will be generated. If this was intentional, please, ignore this "
386 "RESTWarning."
387 << RESTendl;
388 return;
389 }
390
391 string fullPathName = SearchFile((string)fFluxDataFile);
392
393 RESTDebug << "Loading flux table ... " << RESTendl;
394 RESTDebug << "File : " << fullPathName << RESTendl;
395 std::vector<std::vector<Double_t>> fluxData;
396 TRestTools::ReadASCIITable(fullPathName, fluxData, 3);
397
398 RESTDebug << "Table loaded" << RESTendl;
399 TH2F* originalHist = new TH2F("FullTable", "", 100, 0., 1., (Int_t)(20. / fBinSize), 0., 20.);
400 TH2F* continuumHist = new TH2F("ContinuumTable", "", 100, 0., 1., (Int_t)(20. / fBinSize), 0., 20.);
401 TH2F* spectrumHist = new TH2F("LinesTable", "", 100, 0., 1., (Int_t)(20. / fBinSize), 0., 20.);
402
403 Double_t fluxBinSize = TRestTools::GetLowestIncreaseFromTable(fluxData, 1);
404
405 for (const auto& data : fluxData) {
406 Float_t r = 0.005 + data[0];
407 Float_t en = data[1] - 0.005;
408 Float_t flux = data[2] * fluxBinSize; // flux in cm-2 s-1 bin-1
409
410 originalHist->Fill(r, en, (Float_t)flux);
411 continuumHist->Fill(r, en, (Float_t)flux);
412 }
413 RESTDebug << "Histograms filled" << RESTendl;
414
415 Int_t peaks = 0;
416 do {
417 peaks = 0;
418 // We just identify pronounced peaks, not smoothed gaussians!
419 const int smearPoints = (Int_t)(5 / (fBinSize * 100));
420 const int excludePoints = smearPoints / 5;
421 for (const auto& data : fluxData) {
422 Float_t r = 0.005 + data[0];
423 Float_t en = data[1] - 0.005;
424
425 Int_t binR = continuumHist->GetXaxis()->FindBin(r);
426 Int_t binE = continuumHist->GetYaxis()->FindBin(en);
427
428 Double_t avgFlux = 0;
429 Int_t n = 0;
430 for (int e = binE - smearPoints; e <= binE + smearPoints; e++) {
431 if (e < 1 || (e < binE + excludePoints && e > binE - excludePoints)) continue;
432 n++;
433 avgFlux += continuumHist->GetBinContent(binR, e);
434 }
435 avgFlux /= n;
436
437 Float_t targetBinFlux = continuumHist->GetBinContent(binR, binE);
438 Float_t thrFlux = avgFlux + fPeakSigma * TMath::Sqrt(avgFlux);
439 if (targetBinFlux > 0 && targetBinFlux > thrFlux) {
440 continuumHist->SetBinContent(binR, binE, avgFlux);
441 peaks++;
442 }
443 }
444 } while (peaks > 0);
445
446 for (int n = 0; n < originalHist->GetNbinsX(); n++)
447 for (int m = 0; m < originalHist->GetNbinsY(); m++) {
448 Float_t orig = originalHist->GetBinContent(n + 1, m + 1);
449 Float_t cont = continuumHist->GetBinContent(n + 1, m + 1);
450
451 spectrumHist->SetBinContent(n + 1, m + 1, orig - cont);
452 }
453
454 continuumHist->Rebin2D(1, (Int_t)(0.1 / fBinSize)); // cm-2 s-1 (100eV)-1
455 continuumHist->Scale(10); // cm-2 s-1 keV-1
456 // It could be over here if we would use directly a TH2D
457
458 fFluxTable.clear();
459 for (int n = 0; n < continuumHist->GetNbinsX(); n++) {
460 TH1F* hc =
461 (TH1F*)continuumHist->ProjectionY(Form("%s_ContinuumFluxAtRadius%d", GetName(), n), n + 1, n + 1);
462 fFluxTable.push_back(hc);
463 }
464
465 fFluxLines.clear();
466 for (int n = 0; n < spectrumHist->GetNbinsY(); n++) {
467 if (spectrumHist->ProjectionX("", n + 1, n + 1)->Integral() > 0) {
468 Double_t energy = spectrumHist->ProjectionY()->GetBinCenter(n + 1);
469 TH1F* hm = (TH1F*)spectrumHist->ProjectionX(
470 Form("%s_MonochromeFluxAtEnergy%5.3lf", GetName(), energy), n + 1, n + 1);
471 fFluxLines[energy] = hm;
472 }
473 }
474
475 cout << "Number of peaks identified: " << fFluxLines.size() << endl;
476}
477
483 if (fContinuumHist != nullptr) {
484 delete fContinuumHist;
485 fContinuumHist = nullptr;
486 }
487
488 fContinuumHist = new TH1F("ContinuumHist", "", 200, 0, 20);
489 for (const auto& x : fFluxTable) {
490 fContinuumHist->Add(x);
491 }
492
493 fContinuumHist->SetStats(0);
494 fContinuumHist->GetXaxis()->SetTitle("Energy [keV]");
495 fContinuumHist->GetXaxis()->SetTitleSize(0.05);
496 fContinuumHist->GetXaxis()->SetLabelSize(0.05);
497 fContinuumHist->GetYaxis()->SetTitle("Flux [cm-2 s-1 keV-1]");
498 fContinuumHist->GetYaxis()->SetTitleSize(0.05);
499 fContinuumHist->GetYaxis()->SetLabelSize(0.05);
500
501 return fContinuumHist;
502}
503
509 if (fMonoHist != nullptr) {
510 delete fMonoHist;
511 fMonoHist = nullptr;
512 }
513
514 fMonoHist = new TH1F("MonochromaticHist", "", 20000, 0, 20);
515 for (const auto& x : fFluxLines) {
516 fMonoHist->Fill(x.first, x.second->Integral()); // cm-2 s-1 eV-1
517 }
518
519 fMonoHist->SetStats(0);
520 fMonoHist->GetXaxis()->SetTitle("Energy [keV]");
521 fMonoHist->GetXaxis()->SetTitleSize(0.05);
522 fMonoHist->GetXaxis()->SetLabelSize(0.05);
523 fMonoHist->GetYaxis()->SetTitle("Flux [cm-2 s-1 eV-1]");
524 fMonoHist->GetYaxis()->SetTitleSize(0.05);
525 fMonoHist->GetYaxis()->SetLabelSize(0.05);
526
527 return fMonoHist;
528}
529
536 TH1F* hm = GetMonochromaticSpectrum();
537 TH1F* hc = GetContinuumSpectrum();
538
539 if (fTotalHist != nullptr) {
540 delete fTotalHist;
541 fTotalHist = nullptr;
542 }
543
544 fTotalHist = new TH1F("fTotalHist", "", 20000, 0, 20);
545 for (int n = 0; n < hc->GetNbinsX(); n++) {
546 for (int m = 0; m < 100; m++) {
547 fTotalHist->SetBinContent(n * 100 + 1 + m, hc->GetBinContent(n + 1));
548 }
549 }
550
551 for (int n = 0; n < hm->GetNbinsX(); n++)
552 // 1e-2 is the renormalization from 20000 bins to 200 bins
553 fTotalHist->SetBinContent(n + 1, fTotalHist->GetBinContent(n + 1) + 100 * hm->GetBinContent(n + 1));
554
555 fTotalHist->SetStats(0);
556 fTotalHist->GetXaxis()->SetTitle("Energy [keV]");
557 fTotalHist->GetXaxis()->SetTitleSize(0.05);
558 fTotalHist->GetXaxis()->SetLabelSize(0.05);
559 fTotalHist->GetYaxis()->SetTitle("Flux [cm-2 s-1 keV-1]");
560 fTotalHist->GetYaxis()->SetTitleSize(0.05);
561 fTotalHist->GetYaxis()->SetLabelSize(0.05);
562
563 return fTotalHist;
564}
565
570 if (fCanvas != nullptr) {
571 delete fCanvas;
572 fCanvas = nullptr;
573 }
574 fCanvas = new TCanvas("canv", "This is the canvas title", 1200, 500);
575 fCanvas->Draw();
576
577 TPad* pad1 = new TPad("pad1", "This is pad1", 0.01, 0.02, 0.99, 0.97);
578 pad1->Divide(2, 1);
579 pad1->Draw();
580
581 pad1->cd(1);
582 pad1->cd(1)->SetLogy();
583 pad1->cd(1)->SetRightMargin(0.09);
584 pad1->cd(1)->SetLeftMargin(0.15);
585 pad1->cd(1)->SetBottomMargin(0.15);
586
587 TH1F* ht = GetTotalSpectrum();
588 ht->SetLineColor(kBlack);
589 ht->SetFillStyle(4050);
590 ht->SetFillColor(kBlue - 10);
591
592 TH1F* hm = GetMonochromaticSpectrum();
593 hm->SetLineColor(kBlack);
594 hm->Scale(100); // renormalizing per 100eV-1
595
596 ht->Draw("hist");
597 hm->Draw("hist same");
598
599 pad1->cd(2);
600 pad1->cd(2)->SetRightMargin(0.09);
601 pad1->cd(2)->SetLeftMargin(0.15);
602 pad1->cd(2)->SetBottomMargin(0.15);
603
604 ht->Draw("hist");
605 hm->Draw("hist same");
606
607 return fCanvas;
608}
609
615 fFluxLineIntegrals.clear();
617
618 for (const auto& line : fFluxLines) {
619 fTotalMonochromaticFlux += line.second->Integral();
621 }
622
623 for (unsigned int n = 0; n < fFluxLineIntegrals.size(); n++)
625
627 for (unsigned int n = 0; n < fFluxTable.size(); n++) {
628 fTotalContinuumFlux += fFluxTable[n]->Integral() * 0.1; // We integrate in 100eV steps
630 }
631
632 for (unsigned int n = 0; n < fFluxTableIntegrals.size(); n++)
634
636}
637
641Double_t TRestAxionSolarQCDFlux::IntegrateFluxInRange(TVector2 eRange, Double_t mass) {
642 if (eRange.X() == -1 && eRange.Y() == -1) {
644 return GetTotalFlux();
645 }
646
647 Double_t flux = 0;
648 for (const auto& line : fFluxLines)
649 if (line.first > eRange.X() && line.first < eRange.Y()) flux += line.second->Integral();
650
652 for (unsigned int n = 0; n < fFluxTable.size(); n++) {
653 flux += fFluxTable[n]->Integral(fFluxTable[n]->FindFixBin(eRange.X()),
654 fFluxTable[n]->FindFixBin(eRange.Y())) *
655 0.1; // We integrate in 100eV steps
656 }
657
658 return flux;
659}
660
665std::pair<Double_t, Double_t> TRestAxionSolarQCDFlux::GetRandomEnergyAndRadius(TVector2 eRange,
666 Double_t mass) {
667 std::pair<Double_t, Double_t> result = {0, 0};
668 if (!AreTablesLoaded()) return result;
669 Double_t rnd = fRandom->Rndm();
670 if (fTotalMonochromaticFlux == 0 || fRandom->Rndm() > fFluxRatio) {
671 // Continuum
672 for (unsigned int r = 0; r < fFluxTableIntegrals.size(); r++) {
673 if (rnd < fFluxTableIntegrals[r]) {
674 Double_t energy = fFluxTable[r]->GetRandom();
675 if (eRange.X() != -1 && eRange.Y() != -1) {
676 if (energy < eRange.X() || energy > eRange.Y()) return GetRandomEnergyAndRadius(eRange);
677 }
678 Double_t radius = ((Double_t)r + fRandom->Rndm()) * 0.01;
679 std::pair<Double_t, Double_t> p = {energy, radius};
680 return p;
681 }
682 }
683 } else {
684 // Monochromatic
685 int n = 0;
686 for (const auto& line : fFluxLines) {
687 if (rnd < fFluxLineIntegrals[n]) {
688 std::pair<Double_t, Double_t> p = {line.first, line.second->GetRandom()};
689 return p;
690 }
691 n++;
692 }
693 }
694 return result;
695}
696
701 cout << "Continuum solar flux table: " << endl;
702 cout << "--------------------------- " << endl;
703 for (unsigned int n = 0; n < fFluxTable.size(); n++) {
704 for (int m = 0; m < fFluxTable[n]->GetNbinsX(); m++)
705 cout << fFluxTable[n]->GetBinContent(m + 1) << "\t";
706 cout << endl;
707 cout << endl;
708 }
709 cout << endl;
710}
711
716 cout << "Integrated solar flux per solar ring: " << endl;
717 cout << "--------------------------- " << endl;
718 /*
719 for (int n = 0; n < fFluxPerRadius.size(); n++)
720 cout << "n : " << n << " flux : " << fFluxPerRadius[n] << endl;
721 cout << endl;
722 */
723}
724
729 // cout << "Number of monochromatic lines: " << fFluxPerRadius.size() << endl;
730 cout << "+++++++++++++++++++++++++++++++++++" << endl;
731 for (auto const& line : fFluxLines) {
732 cout << "Energy : " << line.first << " keV" << endl;
733 cout << "-----------------" << endl;
734 for (int n = 0; n < line.second->GetNbinsX(); n++)
735 cout << "R : " << line.second->GetBinCenter(n + 1)
736 << " flux : " << line.second->GetBinContent(n + 1) << " cm-2 s-1" << endl;
737 }
738}
739
745
746 if (fFluxDataFile != "")
747 RESTMetadata << " - Solar axion flux datafile (continuum) : " << fFluxDataFile << RESTendl;
748 if (fFluxSptFile != "")
749 RESTMetadata << " - Solar axion flux datafile (monochromatic) : " << fFluxSptFile << RESTendl;
750 RESTMetadata << "-------" << RESTendl;
751 RESTMetadata << " - Total monochromatic flux : " << fTotalMonochromaticFlux << " cm-2 s-1" << RESTendl;
752 RESTMetadata << " - Total continuum flux : " << fTotalContinuumFlux << " cm-2 s-1" << RESTendl;
753 RESTMetadata << "++++++++++++++++++" << RESTendl;
754
759 }
760}
761
767 string rootFilename = TRestTools::GetFileNameRoot(fFluxDataFile);
768
769 string path = REST_USER_PATH + "/export/";
770
771 if (!TRestTools::fileExists(path)) {
772 std::cout << "Creating path: " << path << std::endl;
773 int z = system(("mkdir -p " + path).c_str());
774 if (z != 0) RESTError << "Could not create directory " << path << RESTendl;
775 }
776
777 if (fFluxTable.size() > 0) {
778 std::vector<std::vector<Float_t>> table;
779 for (const auto& x : fFluxTable) {
780 std::vector<Float_t> row;
781 for (int n = 0; n < x->GetNbinsX(); n++) row.push_back(x->GetBinContent(n + 1));
782
783 table.push_back(row);
784 }
785
786 if (!ascii)
787 TRestTools::ExportBinaryTable(path + "/" + rootFilename + ".N200f", table);
788 else
789 TRestTools::ExportASCIITable(path + "/" + rootFilename + ".dat", table);
790 }
791
792 if (fFluxLines.size() > 0) {
793 std::vector<std::vector<Float_t>> table;
794 for (const auto& x : fFluxLines) {
795 std::vector<Float_t> row;
796 row.push_back(x.first);
797 for (int n = 0; n < x.second->GetNbinsX(); n++) row.push_back(x.second->GetBinContent(n + 1));
798
799 table.push_back(row);
800 }
801
803
804 TRestTools::ExportASCIITable(path + "/" + rootFilename + ".spt", table);
805 }
806}
A metadata class to load tabulated solar axion fluxes.
TCanvas * fCanvas
A canvas pointer for drawing.
TRandom3 * fRandom
Random number generator.
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
Bool_t fTablesLoaded
A metadata member to control if this class has been initialized.
A metadata class to load tabulated solar axion fluxes. Mass independent.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarQCDFlux.
std::string fFluxSptFile
The filename containning the solar flux spectra for monochromatic spectrum.
Double_t fPeakSigma
It will be used when loading .flux files to define the threshold for peak identification.
TRestAxionSolarQCDFlux()
Default constructor.
~TRestAxionSolarQCDFlux()
Default destructor.
void ReadFluxFile()
It loads a .flux file. It will split continuum and monochromatic peaks, loading both internal flux ta...
TH1F * GetTotalSpectrum()
It builds a histogram adding the continuum and the monochromatic spectrum component....
TH1F * fMonoHist
A pointer to the monochromatic spectrum histogram.
TH1F * fTotalHist
A pointer to the superposed monochromatic and continuum spectrum histogram.
void LoadContinuumFluxTable()
A helper method to load the data file containning continuum spectra as a function of the solar radius...
void PrintIntegratedRingFlux()
It prints on screen the integrated solar flux per solar ring.
std::vector< Double_t > fFluxTableIntegrals
Accumulative integrated solar flux for each solar ring for continuum spectrum (renormalized to unity)
Double_t GetTotalFlux(Double_t mass=0) override
It returns the total integrated flux at earth in cm-2 s-1.
Double_t IntegrateFluxInRange(TVector2 eRange=TVector2(-1, -1), Double_t mass=0) override
It returns the integrated flux at earth in cm-2 s-1 for the given energy range.
Bool_t LoadTables() override
It defines how to read the solar tables at the inhereted class.
std::vector< TH1F * > fFluxTable
The tabulated solar flux continuum spectra TH1F(200,0,20)keV in cm-2 s-1 keV-1 versus solar radius.
TH1F * GetContinuumSpectrum()
It builds a histogram with the continuum spectrum component. The flux will be expressed in cm-2 s-1 k...
void IntegrateSolarFluxes()
A helper method to initialize the internal class data members with the integrated flux for each solar...
TH1F * fContinuumHist
A pointer to the continuum spectrum histogram.
void LoadMonoChromaticFluxTable()
A helper method to load the data file containning monochromatic spectral lines as a function of the s...
void ExportTables(Bool_t ascii=false) override
It will create files with the continuum and spectral flux components to be used in a later ocasion.
std::vector< Double_t > fFluxLineIntegrals
Accumulative integrated solar flux for each monochromatic energy (renormalized to unity)
TH1F * GetMonochromaticSpectrum()
It builds a histogram with the monochromatic spectrum component. The flux will be expressed in cm-2 s...
std::pair< Double_t, Double_t > GetRandomEnergyAndRadius(TVector2 eRange=TVector2(-1, -1), Double_t mass=0) override
It defines how to generate Monte Carlo energy and radius values to reproduce the flux.
std::map< Double_t, TH1F * > fFluxLines
The tabulated solar flux in cm-2 s-1 for a number of monochromatic energies versus solar radius.
virtual TCanvas * DrawSolarFlux() override
It draws the contents of a .flux file. This method just receives the.
void PrintMonoChromaticFlux()
It prints on screen the spectral lines loaded in memory.
std::string fFluxDataFile
The filename containning the solar flux table with continuum spectrum.
Double_t fFluxRatio
The ratio between monochromatic and total flux.
void PrintContinuumSolarTable()
It prints on screen the table that has been loaded in memory.
Double_t fBinSize
It will be used when loading .flux files to define the input file energy binsize in eV.
Double_t fTotalContinuumFlux
Total solar flux for monochromatic contributions.
Double_t fTotalMonochromaticFlux
Total solar flux for monochromatic contributions.
endl_t RESTendl
Termination flag object for TRestStringOutput.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
std::string SearchFile(std::string filename)
Search files in current directory and directories specified in "searchPath" section.
std::string fConfigFileName
Full name of the rml file.
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
static std::string GetFileNameExtension(const std::string &fullname)
Gets the file extension as the substring found after the latest ".".
Definition: TRestTools.cxx:823
static int ReadASCIITable(std::string fName, std::vector< std::vector< Double_t > > &data, Int_t skipLines=0, std::string separator="\t")
Reads an ASCII file containing a table with values.
Definition: TRestTools.cxx:577
static int ReadBinaryTable(std::string fName, std::vector< std::vector< T > > &data, Int_t columns=-1)
Reads a binary file containing a fixed-columns table with values.
Definition: TRestTools.cxx:253
static std::string GetFileNameRoot(const std::string &fullname)
Gets the filename root as the substring found before the latest ".".
Definition: TRestTools.cxx:834
static T GetLowestIncreaseFromTable(std::vector< std::vector< T > > data, Int_t column)
It returns the lowest increase, different from zero, between the elements of a particular column from...
Definition: TRestTools.cxx:442
static int ExportASCIITable(std::string fname, std::vector< std::vector< T > > &data)
Writes the contents of the vector table given as argument to fname. Allowed types are Int_t,...
Definition: TRestTools.cxx:185
static bool fileExists(const std::string &filename)
Returns true if the file (or directory) with path filename exists.
Definition: TRestTools.cxx:728
static int ExportBinaryTable(std::string fname, std::vector< std::vector< T > > &data)
Writes the contents of the vector table given as argument to fname as a binary file....
Definition: TRestTools.cxx:214
static Bool_t IsBinaryFile(std::string fname)
It identifies if the filename extension follows the formatting ".Nxyzf", where the number of columns ...
Definition: TRestTools.cxx:303
static void TransposeTable(std::vector< std::vector< T > > &data)
It transposes the std::vector<std::vector> table given in the argument. It will transform rows in col...
Definition: TRestTools.cxx:345