CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes | Static Public Attributes | Private Attributes
DQMIO2histo.DQMIO Class Reference

Public Member Functions

def __init__ (self, input_filename, output_filename)
 defined DQMIO types More...
 
def print_index (self)
 
def write_to_file (self, hist_type, index_range, run)
 

Public Attributes

 already_defined
 
 f
 

Static Public Attributes

 types
 

Private Attributes

 _canvas
 
 _filename
 
 _root_file
 

Detailed Description

Class responsible for browsing the content of a DQM file produced
with the DQMIO I/O framework of CMSSW

Definition at line 19 of file DQMIO2histo.py.

Constructor & Destructor Documentation

def DQMIO2histo.DQMIO.__init__ (   self,
  input_filename,
  output_filename 
)

defined DQMIO types

Definition at line 29 of file DQMIO2histo.py.

29  def __init__(self, input_filename, output_filename):
30  self._filename = input_filename
31  self._canvas = None
32  self.f = R.TFile(output_filename, "RECREATE")
33  self.already_defined = {"TProfiles" : False, "TProfile2Ds" : False,
34  "TH2Fs" : False, "TH2Ds" : False}
35 
36  if os.path.exists(self._filename): #we try to open input file if fail
37  self._root_file = R.TFile.Open(self._filename) #-> close script
38  if args.debug:
39  print("## DEBUG ##:")
40  print(" Input: %s\n Output: %s" % (input_filename,
41  output_filename))
42 
43  else:
44  print("File %s does not exists" % self._filename)
45  sys.exit(1)
46 
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:66
def __init__(self, input_filename, output_filename)
defined DQMIO types
Definition: DQMIO2histo.py:29

Member Function Documentation

def DQMIO2histo.DQMIO.print_index (   self)
Loop over the complete index and dump it on the screen.

Definition at line 47 of file DQMIO2histo.py.

References edm.print(), FastTimerService_cff.range, str, DQMIO2histo.DQMIO.types, and DQMIO2histo.DQMIO.write_to_file().

47  def print_index(self):
48  """
49  Loop over the complete index and dump it on the screen.
50  """
51  indices = self._root_file.Get("Indices")
52  if args.debug:
53  print("## DEBUG ##:")
54  print("Run,\tLumi,\tType,\t\tFirstIndex,\tLastIndex")
55  for i in range(indices.GetEntries()):
56  indices.GetEntry(i)
57  print('{0:4d}\t{1:4d}\t{2:4d}({3:s})\t\t{4:4d}\t{5:4d}'.format(
58  indices.Run, indices.Lumi, indices.Type,
59  DQMIO.types[indices.Type], indices.FirstIndex, indices.LastIndex))
60 
61  for i in range(indices.GetEntries()):
62  indices.GetEntry(i)
63  if indices.Type < len(DQMIO.types):
64  self.write_to_file(self.types[indices.Type],
65  [indices.FirstIndex,indices.LastIndex], str(indices.Run))
66 
67  else:
68  print("Unknown histogram type. Type numer: %s" % (indices.Type))
69  self.f.Close()
70 
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:66
def print_index(self)
Definition: DQMIO2histo.py:47
def write_to_file(self, hist_type, index_range, run)
Definition: DQMIO2histo.py:71
#define str(s)
def DQMIO2histo.DQMIO.write_to_file (   self,
  hist_type,
  index_range,
  run 
)
Method looping over entries for specified histogram type and 
writing to FullName path to output ROOT File

Definition at line 71 of file DQMIO2histo.py.

References DQMIO2histo.DQMIO.already_defined, join(), edm.print(), FastTimerService_cff.range, and str.

Referenced by DQMIO2histo.DQMIO.print_index().

71  def write_to_file(self, hist_type, index_range, run):
72  """
73  Method looping over entries for specified histogram type and
74  writing to FullName path to output ROOT File
75  """
76  print("Working on: %s indexes: %s..%s" % (hist_type ,index_range[0],
77  index_range[1]))
78  t_tree = self._root_file.Get(hist_type)
79  __run_dir = "Run %s" % (run)
80  ###we set Branch for the needed type
81  if hist_type == "TProfiles":
82  if not self.already_defined["TProfiles"]:
83  R.gROOT.ProcessLine("TProfile* _tprof;")
84  self.already_defined["TProfiles"] = True
85  t_tree.SetBranchAddress("Value", R._tprof)
86  t_tree.GetEntry(index_range[0])
87  elif hist_type == "TProfile2Ds":
88  if not self.already_defined["TProfile2Ds"]:
89  R.gROOT.ProcessLine("TProfile2D* _tprof2d;")
90  self.already_defined["TProfile2Ds"] = True
91  t_tree.SetBranchAddress("Value", R._tprof2d)
92  t_tree.GetEntry(index_range[0])
93  elif hist_type == "TH2Fs":
94  if not self.already_defined["TH2Fs"]:
95  R.gROOT.ProcessLine("TH2F* _th2f;")
96  self.already_defined["TH2Fs"] = True
97  t_tree.SetBranchAddress("Value", R._th2f)
98  t_tree.GetEntry(index_range[0])
99  elif hist_type == "TH2Ds":
100  if not self.already_defined["TH2Ds"]:
101  R.gROOT.ProcessLine("TH2D* _th2d;")
102  self.already_defined["TH2Ds"] = True
103  t_tree.SetBranchAddress("Value", R._th2d)
104  t_tree.GetEntry(index_range[0])
105 
106  for i in range(0,t_tree.GetEntries()+1): ##iterate on entries for specified type
107  if i >= index_range[0] and i <= index_range[1]: ##if entries as in range:
108  t_tree.GetEntry(i)
109  name = str(t_tree.FullName)
110  # print " %s: %s" % (i, name)
111  file_path = name.split("/")[:-1] ## same run/lumi histograms
112  __directory = "%s/%s" % (os.path.join("DQMData", __run_dir),
113  "/".join(file_path))
114  directory_ret = self.f.GetDirectory(__directory)
115  if not directory_ret:
116  self.f.mkdir(os.path.join(__directory))
117  self.f.cd(os.path.join(__directory))
118  if hist_type == "Strings":
119  construct_str = '<%s>s=%s</%s>' % (name.split("/")[-1:][0],
120  t_tree.Value, name.split("/")[-1:][0])
121  tmp_str = R.TObjString(construct_str)
122  tmp_str.Write()
123  elif hist_type == "Ints":
124  construct_str = '<%s>i=%s</%s>' % (name.split("/")[-1:][0],
125  t_tree.Value, name.split("/")[-1:][0])
126  tmp_str = R.TObjString(construct_str)
127  tmp_str.Write()
128  elif hist_type == "Floats":
129  construct_str = '<%s>f=%s</%s>' % (name.split("/")[-1:][0],
130  t_tree.Value, name.split("/")[-1:][0])
131  tmp_str = R.TObjString(construct_str)
132  tmp_str.Write()
133  else:
134  if hist_type in ["TProfiles", "TProfile2Ds", "TH2Fs", "TH2Ds"]:
135  if hist_type == "TProfiles": #if type is specific we write it.
136  R._tprof.Write()
137  elif hist_type == "TProfile2Ds":
138  R._tprof2d.Write()
139  elif hist_type == "TH2Fs":
140  R._th2f.Write()
141  elif hist_type == "TH2Ds":
142  R._th2d.Write()
143  else: #else we wirte Leafs Value which is a histogram
144  t_tree.Value.Write()
145 
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:66
static std::string join(char **cmd)
Definition: RemoteFile.cc:17
def write_to_file(self, hist_type, index_range, run)
Definition: DQMIO2histo.py:71
#define str(s)

Member Data Documentation

DQMIO2histo.DQMIO._canvas
private

Definition at line 31 of file DQMIO2histo.py.

DQMIO2histo.DQMIO._filename
private
DQMIO2histo.DQMIO._root_file
private

Definition at line 37 of file DQMIO2histo.py.

DQMIO2histo.DQMIO.already_defined

Definition at line 33 of file DQMIO2histo.py.

Referenced by DQMIO2histo.DQMIO.write_to_file().

DQMIO2histo.DQMIO.f
DQMIO2histo.DQMIO.types
static

Definition at line 24 of file DQMIO2histo.py.

Referenced by DQMIO2histo.DQMIO.print_index().