hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
SvtBlFitHistoProcessor.cxx
Go to the documentation of this file.
2#include <FlatTupleMaker.h>
3#include <string>
4#include <algorithm>
5#include <cstdlib>
7 : Processor(name, process) {
8 }
9
12
14
15 std::cout << "[SvtBlFitHistoProcessor] Configuring" << std::endl;
16 try
17 {
18 histCfgFilename_ = parameters.getString("histCfg");
19 rawhitsHistCfgFilename_ = parameters.getString("rawhitsHistCfg");
20 thresholdsFileIn_ = parameters.getString("thresholdsFileIn");
21 layer_ = parameters.getString("layer");
22 rebin_ = parameters.getInteger("rebin");
23 minStats_ = parameters.getInteger("minStats");
24 deadRMS_ = parameters.getInteger("deadRMS");
25 debug_ = parameters.getInteger("debug");
26 year_ = parameters.getInteger("year");
27 }
28 catch (std::runtime_error& error)
29 {
30 std::cout << error.what() << std::endl;
31 }
32}
33
34void SvtBlFitHistoProcessor::initialize(std::string inFilename, std::string outFilename) {
35
36 std::cout << "initializing SvtBlFitHistoProcessor" << std::endl;
37 //InFile containing 2D histograms from the SvtCondAnaProcessor
38 inF_ = new TFile(inFilename.c_str());
39 outF_ = new TFile(outFilename.c_str(),"RECREATE");
40
41 //init tuple to store fit values
42 flat_tuple_ = new FlatTupleMaker(outFilename.c_str(), "gaus_fit");
43
44 //Initialize fit histos
47 std::cout << "[BlFitHistos] Loading 2D Histos" << std::endl;
50
51
52 //Setup flat tuple branches
53 //Name of each hybrid
54 flat_tuple_->addString("halfmodule_hh");
55 //number of entries in channel
56 flat_tuple_->addString("n_entries");
57 //minimum required entries in channel to attempt fit
58 flat_tuple_->addVariable("minStats");
59 //channel rebinning setting
60 flat_tuple_->addVariable("rebin");
61 //channel number
62 flat_tuple_->addVariable("channel");
63 //svt_id
64 flat_tuple_->addVariable("svt_id");
65 //if 1.0, channel stats too low to perform fit
66 flat_tuple_->addVariable("lowStats");
67 //channel rms. Use to flag dead channels
69 //if rms < threshold, channel is dead
70 flat_tuple_->addVariable("dead");
71 //if channel fit is grossly poor
72 flat_tuple_->addVariable("badfit");
73 //if fit mean + N*sigma > xmax, flag as lowdaq
74 flat_tuple_->addVariable("lowdaq");
75 //if fit mean + N*sigma > xmax, flag as lowdaq
76 flat_tuple_->addVariable("suplowDaq");
77 //Threshold value setting on apv channel during run
78 flat_tuple_->addVariable("threshold");
79 flat_tuple_->addVariable("minthreshold");
80
81 //Fit values
82 flat_tuple_->addVariable("BlFitMean");
83 flat_tuple_->addVariable("BlFitNorm");
84 flat_tuple_->addVariable("BlFitSigma");
85 flat_tuple_->addVariable("BlFitRangeLower");
86 flat_tuple_->addVariable("BlFitRangeUpper");
87 flat_tuple_->addVariable("BlFitChi2");
88 flat_tuple_->addVariable("BlFitNdf");
89}
90
91
93 if(debug_)
94 std::cout << "RUNNING PROCESS ON 2d HISTOS" << std::endl;
95 std::map<std::string,TH2F*> histos2d = fitHistos_->get2dHistos();
96 std::map<std::string, TH2F*>::iterator it;
97 for(it = histos2d.begin(); it != histos2d.end(); it++)
98 std::cout << "2d Histo Loaded: " << it->first << std::endl;
99 if(debug_)
100 std::cout << "FIT 2d HISTO Baselines" << std::endl;
102 return true;
103}
104
106
107 std::cout << "finalizing SvtBlFitHistoProcessor" << std::endl;
109 outF_->Close();
110}
111
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
description
Definition BlFitHistos.h:30
void setDebug(bool value)
Set debug.
Definition BlFitHistos.h:91
std::map< std::string, TH2F * > get2dHistos()
description
Definition BlFitHistos.h:50
void getHistosFromFile(TFile *inFile, std::string layer="")
Get the Histos From File object.
void fit2DHistoChannelBaselines(std::map< std::string, TH2F * > histos2d, int rebin_, int minStats_, int deadRMS_, std::string thresholdsFileIn_, FlatTupleMaker *flat_tuple_)
description
description
void addVariable(std::string variable_name)
description
void close()
description
void addString(std::string variable_name)
description
virtual void loadHistoConfig(const std::string histoConfigFile)
load histogram config
description
Base class for all event processing components.
Definition Processor.h:34
TFile * outF_
Definition Processor.h:125
Insert description here. more details.
virtual void configure(const ParameterSet &parameters)
description  *  *
virtual void finalize()
description  * 
SvtBlFitHistoProcessor(const std::string &name, Process &process)
Constructor  *  *.
std::string layer_
description
std::string thresholdsFileIn_
description
BlFitHistos * fitHistos_
description
std::string histCfgFilename_
description
virtual void initialize(std::string inFilename, std::string outFilename)
description  *  *
FlatTupleMaker * flat_tuple_
description
std::string rawhitsHistCfgFilename_
description
virtual bool process()
description  *  *