Jpp  19.1.0
the software that should make you happy
Functions
JMultiplicityK40.cc File Reference

Example program to calculate multiples rate. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <map>
#include <cmath>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TRandom3.h"
#include "JTools/JRange.hh"
#include "JTools/JWeight.hh"
#include "JTools/JQuantile.hh"
#include "JPhysics/KM3NeT.hh"
#include "JPhysics/KM3NeT2D.hh"
#include "JPhysics/JConstants.hh"
#include "JGeometry3D/JPosition3D.hh"
#include "JGeometry3D/JDirection3D.hh"
#include "JMath/JMathToolkit.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDetectorSupportkit.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Detailed Description

Example program to calculate multiples rate.

Author
mdejong

Definition in file JMultiplicityK40.cc.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 682 of file JMultiplicityK40.cc.

683 {
684  using namespace std;
685  using namespace JPP;
686 
687  typedef long long int counter_type;
688 
689  string outputFile;
690  double bequerel;
691  JRange_t D_m;
692  counter_type numberOfEvents;
693  double QE;
694  UInt_t seed;
695  int generator;
696  double focus;
697  int Aefftype;
698  bool exclusive;
699  int debug;
700 
701  try {
702 
703  JParser<> zap("Example program to calculate multiples rate.");
704 
705  zap['o'] = make_field(outputFile) = "k40.root";
706  zap['b'] = make_field(bequerel) = 13.75e3; // [m^-3 s^-1]
707  zap['n'] = make_field(numberOfEvents) = 1e7;
708  zap['D'] = make_field(D_m) = JRange_t(0.216, 10);
709  zap['Q'] = make_field(QE) = 1.0;
710  zap['S'] = make_field(seed) = 0;
711  zap['G'] = make_field(generator) = 0, +2, -2;
712  zap['F'] = make_field(focus) = 1.0;
713  zap['A'] = make_field(Aefftype) = 1, 2, 3;
714  zap['U'] = make_field(exclusive);
715  zap['a'] = make_field(a) = 0.0;
716  zap['d'] = make_field(debug) = 2;
717 
718  zap(argc, argv);
719  }
720  catch(const exception &error) {
721  FATAL(error.what() << endl);
722  }
723 
724  gRandom->SetSeed(seed);
725 
726  using namespace NAMESPACE;
727 
728 
729  const int id = 1;
730  const JModule module = getModule<JKM3NeT_t>(id);
731 
732  DEBUG(module << endl);
733 
734  bequerel /= focus; // correct decay rate for focussing of photons
735 
736  const double R_m = 17.0 * 2.54 * 0.5e-2; // radius 17'' optical module [m]
737  const double A = PI * R_m * R_m; // cross section optical module [m^2]
738 
739  const double wmin = 280.0; // minimal wavelength [nm]
740  const double wmax = 700.0; // maximal wavelength [nm]
741  double ng = 41.0; // average number of photons per decay in given wavelength range
742 
743  const double WAVELENGTH_EXPANSION = (wmax-wmin) / (wmin*wmax) * (300.0*600.0)/(600.0-300.0);
744 
745  JGenerator* enigma = NULL; // distance generation
746 
747  switch (generator) {
748 
749  case +2: enigma = new JEnigma<+2>(D_m); break;
750  case 0: enigma = new JEnigma< 0>(D_m); break;
751  case -2: enigma = new JEnigma<-2>(D_m); break;
752 
753  default: FATAL("No generator type " << generator << endl);
754  }
755 
756  const double vmin = 1.0 / wmax; // wavelength generation
757  const double vmax = 1.0 / wmin; // as lambda^-2
758 
759  double QEmax = 0.0;
760 
761  for (double w = wmin; w <= wmax; w += 1.0) {
762  if (getQE(w) > QEmax) {
763  QEmax = getQE(w);
764  }
765  }
766 
767  NOTICE("Maximal QE " << FIXED(5,3) << QEmax << endl);
768  NOTICE("Wavelength expansion " << FIXED(5,3) << WAVELENGTH_EXPANSION << endl);
769  NOTICE("Number of photons per decay " << FIXED(5,2) << ng << endl);
770 
771  typedef JManager<int, TH1D> JManager_t;
772 
773  JManager_t H1(new TH1D("M[%]", NULL, 100, D_m.getLowerLimit(), D_m.getUpperLimit()));
774 
775  H1->Sumw2();
776 
777  TH1D pmt("pmt", NULL, 1000, -1.0, +1.0);
778 
779  for (Int_t i = 1; i != pmt.GetNbinsX(); ++i) {
780 
781  const double dot = pmt.GetBinCenter(i);
782 
783  double y = 0.0;
784 
785  switch (Aefftype) {
786 
787  case 1:
788  y = getAngularAcceptance(dot);
789  break;
790 
791  case 3:
792  y = get_angular_acceptance(dot);
793  break;
794  }
795 
796  pmt.SetBinContent(i, y);
797  }
798 
799 
802 
803 
804  vector<double> pi(module.size());
805  vector<int> buffer;
806 
807  for (counter_type event_count = 0; event_count != numberOfEvents; ++event_count) {
808 
809  if (event_count%10000 == 0) {
810  STATUS("event: " << setw(10) << event_count << '\r'); DEBUG(endl);
811  }
812 
813  const JResult& result = enigma->next();
814 
815  const double D = result.D;
816  const double V = result.V;
817 
818  // check first photons on cross section of optical module at its face
819  // it is assumed that the light from the K40 decay is purely isotropic
820 
821  double W = A / (4*PI*(D-R_m)*(D-R_m));
822 
823  if (W > 0.5) {
824  W = 0.5;
825  }
826 
827  //
828  double x = gRandom->Rndm(); // decay mode
829  double y = 0.0; // number of photons
830 
831  if ((x -= k40_beta_decay .getBranchingRatio()) <= 0.0)
832  y = k40_beta_decay (gRandom->Rndm());
833  else if ((x -= k40_electron_capture.getBranchingRatio()) <= 0.0)
834  y = k40_electron_capture(gRandom->Rndm());
835 
836  const int N = gRandom->Poisson(y * WAVELENGTH_EXPANSION * QE * W * QEmax * focus);
837  //
838  /*
839  const int N = gRandom->Poisson(ng * WAVELENGTH_EXPANSION * QE * W * QEmax * focus);
840  */
841  if (N >= 2) {
842 
843  // generate vertex
844 
845  const double ct = gRandom->Uniform(-1.0, +1.0);
846  const double phi = gRandom->Uniform(-PI, +PI);
847 
848  const double st = sqrt((1.0 - ct) * (1.0 + ct));
849 
850  const JPosition3D vertex(D * st * cos(phi),
851  D * st * sin(phi),
852  D * ct);
853 
854  buffer.clear();
855 
856  for (int i = 0; i != N; ++i) {
857 
858  // generate wavelength of photon
859 
860  const double v = gRandom->Uniform(vmin, vmax);
861  const double w = 1.0 / v;
862 
863  const double l_abs = getAbsorptionLength(w);
864 
865  double P = 0.0;
866 
867  for (size_t pmt = 0; pmt != module.size(); ++pmt) {
868 
869  JPosition3D pos(module[pmt].getPosition());
870 
871  pos -= vertex;
872 
873  const double d = pos.getLength();
874 
875  pos /= d; // direction of photon
876 
877  if (d < D - R_m) {
878  ERROR("Distance " << d << " < " << D << endl);
879  }
880 
881  const double dot = getDot(pos, module[pmt].getDirection());
882 
883  const double U = 0.5 * (1.0 - d / sqrt(d*d + getPhotocathodeArea() / PI));
884 
885  double p = 0.0;
886 
887  switch (Aefftype) {
888 
889  case 1:
890  p = getAngularAcceptance(dot) * getQE(w);
891  break;
892 
893  case 2:
895  break;
896 
897  case 3:
898  p = get_angular_acceptance(dot) * getQE(w);
899  break;
900  }
901 
902  P += pi[pmt] = U * p * exp(-d/l_abs);
903  }
904 
905  if (P > W) {
906  ERROR("Probability " << P << " > " << W << endl);
907  }
908 
909  if (W * QEmax * gRandom->Rndm() < P) {
910 
911  int pmt = 0;
912  double y = gRandom->Uniform(P);
913 
914  for (vector<double>::const_iterator i = pi.begin(); i != pi.end() && (y -= *i) > 0.0; ++i, ++pmt) {}
915 
916  buffer.push_back(pmt);
917  }
918  }
919 
920  if (!buffer.empty()) {
921 
922  int M = buffer.size();
923 
924  if (exclusive) {
925 
926  sort(buffer.begin(), buffer.end());
927 
928  M = distance(buffer.begin(), unique(buffer.begin(), buffer.end()));
929  }
930 
931  h1[M].put(V);
932  H1[M]->Fill(D, V);
933 
934  for (int i = 2; i <= M; ++i) {
935  P2[i].put((double) (buffer.size() - M) / (double) M, V);
936  }
937 
938  }
939  }
940  }
941  STATUS(endl);
942 
943  for (JManager_t::iterator i = H1.begin(); i != H1.end(); ++i) {
944  i->second->Scale(bequerel / (double) numberOfEvents);
945  }
946 
947  for (size_t M = 2; M != 7; ++M) {
948  cout << "Rate[" << M << "] = "
949  << FIXED(8,3) << bequerel * h1[M].getTotal() / (double) numberOfEvents
950  << " +/- "
951  << FIXED(7,3) << bequerel * h1[M].getError() / (double) numberOfEvents
952  << " Hz" << endl;
953  }
954 
955  for (size_t M = 2; M != 7; ++M) {
956  if (P2[M].getCount() != 0) {
957  cout << "P2[" << M << "] = " << P2[M].getMean() << endl;
958  }
959  }
960 
961  TFile out(outputFile.c_str(), "recreate");
962 
963  out << H1 << pmt;
964 
965  out.Write();
966  out.Close();
967 }
string outputFile
double getAngularAcceptance(const double x)
Angular acceptence of PMT.
Definition: JDrawLED.cc:68
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define STATUS(A)
Definition: JMessage.hh:63
#define ERROR(A)
Definition: JMessage.hh:66
#define NOTICE(A)
Definition: JMessage.hh:64
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Data structure for a composite optical module.
Definition: JModule.hh:75
Data structure for position in three dimensions.
Definition: JPosition3D.hh:38
Utility class to parse command line options.
Definition: JParser.hh:1698
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
double getPhotocathodeArea()
Get photo-cathode area of PMT.
Definition: Antares.hh:51
double getAbsorptionLength(const double lambda)
Get absorption length.
Definition: Antares.hh:63
const double a
Definition: JQuadrature.cc:42
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
double getDot(const JNeutrinoDirection &first, const JNeutrinoDirection &second)
Dot product.
Definition: JAstronomy.hh:676
double getQE(const double R, const double mu)
Get QE for given ratio of hit probabilities and expectation value of the number of photo-electrons.
size_t getCount(const array_type< T > &buffer, const JCompare_t &compare)
Count number of unique values.
Definition: JVectorize.hh:261
static const double PI
Mathematical constants.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Long64_t counter_type
Type definition for counter.
data_type w[N+1][M+1]
Definition: JPolint.hh:867
data_type v[N+1][M+1]
Definition: JPolint.hh:866
static const JPhotocathodeArea2D getPhotocathodeArea2D
Function object for effective photo-cathode area of PMT.
Definition: KM3NeT2D.hh:5235
Definition: JSTDTypes.hh:14
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:448
Description of Monte Carlo event generation applications.
Definition: JHead.hh:469