#include #include #include #include "TROOT.h" #include "TFile.h" #include "TProfile.h" #include "TProfile2D.h" #include "JDetector/JPMTSimulator.hh" #include "JDetector/JPMTAnalogueSignalProcessor.hh" #include "JDetector/JPMTDefaultSimulator.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" /** * \file * * Example program to histogram PMT simulator features. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; string outputFile; int numberOfHits; JPMTParameters parameters; int debug; try { JProperties properties = parameters.getProperties(); JParser<> zap("Example program to histogram PMT simulator features."); zap['o'] = make_field(outputFile) = "pmt.root"; zap['n'] = make_field(numberOfHits) = 100; zap['P'] = make_field(properties) = JPARSER::initialised(); zap['d'] = make_field(debug) = 0; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } if (debug >= JEEP::debug_t) { cout << "PMT parameters:" << endl; cout << parameters.getProperties(JEquationParameters("=", "\n", "", "")) << endl; } const JPMTAnalogueSignalProcessor cpu(parameters); const JPMTIdentifier pmt(1,0); const JCalibration calibration; const JStatus status; const JPMTDefaultSimulator simulator(parameters, pmt); JPMTData input; JPMTData output; TFile out(outputFile.c_str(), "recreate"); TProfile h0("tot", NULL, 1000, parameters.threshold, 1000.0); TProfile h1("hit", NULL, 1000, parameters.threshold, 1000.0); TProfile2D h2("2D", NULL, 50, 0.0, 250.0, 100, 0.0, 1000.0); for (int i = 1; i <= h1.GetNbinsX(); ++i) { const double t_ns = 100.0; const double x = h1.GetBinCenter(i); const int npe = (int) x; for (int i = 0; i != numberOfHits; ++i) { input .clear(); output.clear(); input.push_back(JPMTSignal(t_ns, npe)); simulator.processHits(pmt, calibration, status, input, output); double tot = 0.0; for (JPMTData::const_iterator hit = output.begin(); hit != output.end(); ++hit) { tot += hit->tot_ns; } h0.Fill(x, tot); h1.Fill(x, (double) output.size()); } } for (int i = 1; i <= h2.GetXaxis()->GetNbins(); ++i) { for (int j = 1; j <= h2.GetYaxis()->GetNbins(); ++j) { const double t_ns = (int) h2.GetXaxis()->GetBinCenter(i); const int npe = (int) h2.GetYaxis()->GetBinCenter(j); for (int i = 0; i != numberOfHits; ++i) { input .clear(); output.clear(); input.push_back(JPMTSignal(0.0, npe)); input.push_back(JPMTSignal(t_ns, 1)); simulator.processHits(pmt, calibration, status, input, output); if (!output.empty()) { h2.Fill(t_ns, (double) npe, output.begin()->tot_ns); } } } } out.Write(); out.Close(); }