MARLEY (Model of Argon Reaction Low Energy Yields)  v1.2.0
A Monte Carlo event generator for tens-of-MeV neutrino interactions
marley_root.hh
1 //
5 // This file is part of MARLEY (Model of Argon Reaction Low Energy Yields)
6 //
7 // MARLEY is free software: you can redistribute it and/or modify it under the
8 // terms of version 3 of the GNU General Public License as published by the
9 // Free Software Foundation.
10 //
11 // For the full text of the license please see COPYING or
12 // visit http://opensource.org/licenses/GPL-3.0
13 //
14 // Please respect the MCnet academic usage guidelines. See GUIDELINES
15 // or visit https://www.montecarlonet.org/GUIDELINES for details.
16 
17 #pragma once
18 
19 // standard library includes
20 #include <memory>
21 
22 // ROOT includes
23 #include "TFile.h"
24 #include "TGraph.h"
25 #include "TH1.h"
26 
27 // MARLEY includes
28 #include "marley/Error.hh"
29 #include "marley/NeutrinoSource.hh"
30 
31 namespace marley_root {
32 
40  template<typename T> T* get_root_object(const std::string& tfile_name,
41  const std::string& namecycle);
42 
43  std::unique_ptr<marley::GridNeutrinoSource>
44  make_root_neutrino_source(int pdg_code, const TH1* th1);
45 
46  std::unique_ptr<marley::GridNeutrinoSource>
47  make_root_neutrino_source(int pdg_code, const TGraph* tgraph);
48 }
49 
50 template<typename T> T*
51  marley_root::get_root_object(const std::string& tfile_name,
52  const std::string& namecycle)
53 {
54  // Attempt to open the TFile and complain if something goes wrong.
55  // Use a std::unique_ptr so that the TFile object will be auto-deleted
56  // (and therefore closed) when this function terminates.
57  std::unique_ptr<TFile> file(TFile::Open(tfile_name.c_str(),
58  "read"));
59 
60  if (file) {
61  // Attempt to read in the ROOT object from the file using the given
62  // namecycle. Return pointer to it (or nullptr if something went wrong)
63  T* obj = dynamic_cast<T*>(file->Get(namecycle.c_str()));
64  // Force the TFile to disown the object if it inherits from TH1
65  // (otherwise, ROOT will auto-delete it when the TFile object is
66  // deleted).
67  // TODO: decide what to do if the type is TTree (disown too? throw
68  // exception?) since TTrees also behave this way
69  // TODO: add other checks if you discover more ROOT classes that force
70  // ownership
71  TH1* th1_test = dynamic_cast<TH1*>(obj);
72  if (th1_test) th1_test->SetDirectory(nullptr);
73  // Return a std::unique_ptr to the object
74  return dynamic_cast<T*>(obj);
75  }
76 
77  // couldn't open ROOT file
78  else throw marley::Error(std::string("Failed to open") + " ROOT file '"
79  + tfile_name + "'");
80 
81  return nullptr;
82 }
Base class for all exceptions thrown by MARLEY functions.
Definition: Error.hh:32