#include #include #include #include #include #include "TROOT.h" #include "TFile.h" #include "TH2D.h" #include "JTools/JFunction1D_t.hh" #include "JTools/JFunctionalMap_t.hh" #include "JTools/JAbstractHistogram.hh" #include "JPhysics/JPDFTransformer.hh" #include "JPhysics/JPDFTable.hh" #include "JPhysics/JPDFTypes.hh" #include "JPhysics/JNPETable.hh" #include "JPhysics/JGeanz.hh" #include "JGeometry3D/JAngle3D.hh" #include "Jeep/JPrint.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" /** * \file * * Program to plot PDF of Cherenkov light from EM-shower using interpolation tables. * \author lquinn */ int main(int argc, char **argv) { using namespace std; using namespace JPP; typedef JAbstractHistogram JHistogram_t; vector inputFile; string outputFile; JAngle3D dir; JHistogram_t x; JHistogram_t y; double E; int numberOfPoints; int debug; try { JParser<> zap("Program to plot PDF of Cherenkov light from EM-shower using interpolation tables."); zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = ""; zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]"); zap['x'] = make_field(x, "histogram x-binning") = JHistogram_t(); zap['y'] = make_field(y, "histogram y-binning") = JHistogram_t(); zap['E'] = make_field(E, "Energy [GeV]") = 0.0; zap['N'] = make_field(numberOfPoints) = 10; zap['d'] = make_field(debug) = 0; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } vector Z; if (E > 0.0) { for (int i = 0; i != numberOfPoints; ++i) { Z.push_back(geanz.getLength(E, (i + 0.5) / (double) numberOfPoints)); } } typedef JPolint0Function1D_t JFunction1D_t; typedef JMapList > > > JMapList_t; typedef JPDFTable JPDF_t; typedef JNPETable JNPE_t; JDistance::precision = 1.0e-10; const int N = inputFile.size(); JPDF_t pdf[N]; JNPE_t npe[N]; try { for (int i = 0; i != N; ++i) { NOTICE("loading input from file " << inputFile[i] << "... " << flush); pdf[i].load(inputFile[i].c_str()); pdf[i].setExceptionHandler(new JFunction1D_t::JDefaultResult(JMATH::zero)); npe[i] = JNPE_t(pdf[i]); NOTICE("OK" << endl); } } catch(const JException& error) { FATAL(error.what() << endl); } TFile out(outputFile.c_str(), "recreate"); if (!x.is_valid()) { x = JHistogram_t(150, 0.0, 150.0); } if (!y.is_valid()) { y = JHistogram_t(200, -1.0, +1.0); } TH2D h0("h0", "PDF Projection; D [m]; cos #theta_{0}", x.getNumberOfBins(), x.getLowerLimit(), x.getUpperLimit(), y.getNumberOfBins(), y.getLowerLimit(), y.getUpperLimit()); for (int ix = 1; ix <= h0.GetNbinsX(); ++ix) { for (int iy = 1; iy <= h0.GetNbinsY(); ++iy) { const double cd = h0.GetYaxis()->GetBinCenter(iy); const double D = h0.GetXaxis()->GetBinCenter(ix); double Y = 0.0; if (!Z.empty()) { const double W = 1.0 / (double) Z.size(); for (vector::const_iterator z = Z.begin(); z != Z.end(); ++z) { const double __D = sqrt(D*D - 2.0*(D*cd)*(*z) + (*z)*(*z)); const double __cd = (D * cd - (*z)) / __D; for (int i = 0; i != N; ++i) { try { Y += W * npe[i](__D, __cd, dir.getTheta(), dir.getPhi()); } catch(const exception& error) {} } } } else { for (int i = 0; i != N; ++i) { Y += npe[i](D, cd, dir.getTheta(), dir.getPhi()); } } h0.SetBinContent(ix, iy, Y); } } out.Write(); out.Close(); }