ATLAS Offline Software
Loading...
Searching...
No Matches
ToolsSplit.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
5//
6// Distributed under the Boost Software License, Version 1.0.
7// (See accompanying file LICENSE_1_0.txt or copy at
8// http://www.boost.org/LICENSE_1_0.txt)
9
10// Please feel free to contact me (krumnack@iastate.edu) for bug
11// reports, feature suggestions, praise and complaints.
12
13
14//
15// includes
16//
17
19
27#include <TFile.h>
28#include <TTree.h>
29#include <memory>
30#include <sstream>
31
32//
33// method implementations
34//
35
36namespace SH
37{
39 {
40 for (SampleHandler::iterator sample = sh.begin(),
41 end = sh.end(); sample != end; ++ sample)
42 scanNEvents (**sample);
43 }
44
45
46
47 void scanNEvents (Sample& sample)
48 {
49 SampleLocal *const mysample = dynamic_cast<SampleLocal*>(&sample);
50 if (!mysample)
51 RCU_THROW_MSG ("sample not of type SampleLocal");
52
53 const std::string tree_name = sample.meta()->castString (MetaFields::treeName, MetaFields::treeName_default);
54 if (tree_name.empty())
55 RCU_THROW_MSG ("sample doesn't contain a tree name");
56
57 Long64_t tot_entries = 0;
58 std::vector<Long64_t> entries;
59 const auto fileList = sample.makeFileList ();
60 for (const std::string& file_name : fileList)
61 {
62 std::unique_ptr<TFile> file (TFile::Open (file_name.c_str(), "READ"));
63 if (!file.get())
64 RCU_THROW_MSG ("failed to open file " + file_name);
65 Long64_t treeEntries = 0;
66 TTree *const tree = dynamic_cast<TTree*>(file->Get (tree_name.c_str()));
67 if (tree != 0)
68 treeEntries = tree->GetEntries ();
69 entries.push_back (treeEntries);
70 tot_entries += treeEntries;
71 }
72 RCU_ASSERT (entries.size() == fileList.size());
73 sample.meta()->addReplace (new MetaVector<Long64_t>(MetaFields::numEventsPerFile, entries));
74 sample.meta()->setDouble (MetaFields::numEvents, tot_entries);
75 }
76
77
78
79 SampleHandler splitSample (Sample& sample, const Long64_t nevt)
80 {
81 if (!dynamic_cast<SampleLocal*>(&sample))
82 RCU_THROW_MSG ("sample not of type SampleLocal");
83
84 TObject *meta = sample.meta()->get (MetaFields::numEventsPerFile);
85 if (!meta)
86 {
87 RCU_WARN_MSG ("sample " + sample.name() + " lacks nc_nevtfile, running scanNEvents, please save sample");
88 scanNEvents (sample);
89 meta = sample.meta()->get (MetaFields::numEventsPerFile);
90 }
91 RCU_ASSERT (meta != 0);
92 MetaVector<Long64_t> *const nentries
93 = dynamic_cast<MetaVector<Long64_t> *>(meta);
94 if (nentries == 0)
95 RCU_THROW_MSG ("nc_nevtfile is of the wrong type");
96 if (nentries->value.size() != sample.numFiles())
97 RCU_THROW_MSG ("nc_nevtfile has the wrong number of entries");
98
100 std::unique_ptr<SampleLocal> res;
101 Long64_t num = 0;
102 const std::string meta_tree = sample.meta()->castString (MetaFields::treeName, MetaFields::treeName_default);
103 const double meta_xs = sample.meta()->castDouble (MetaFields::crossSection, 0);
104 for (unsigned file = 0, end = nentries->value.size(); file != end; ++ file)
105 {
106 if (num > 0 && num + nentries->value[file] > nevt)
107 {
108 result.add (res.release());
109 num = 0;
110 }
111 if (res.get() == 0)
112 {
113 std::ostringstream name;
114 name << sample.name() << "_" << result.size();
115 res.reset (new SampleLocal (name.str()));
116 res->tags (sample.tags());
117 res->meta()->fetch (*sample.meta());
118 if (!meta_tree.empty())
119 res->meta()->setString (MetaFields::treeName, meta_tree);
120 if (meta_xs != 0)
121 res->meta()->setDouble (MetaFields::crossSection, meta_xs);
122 }
123 res->add (sample.fileName (file));
124 num += nentries->value[file];
125 }
126 if (num > 0)
127 result.add (res.release());
128 return result;
129 }
130}
#define RCU_ASSERT(x)
Definition Assert.h:222
std::pair< std::vector< unsigned int >, bool > res
#define RCU_THROW_MSG(message)
Definition PrintMsg.h:58
#define RCU_WARN_MSG(message)
Definition PrintMsg.h:52
This class defines a templatized version of the meta-data in vector form.
Definition MetaVector.h:28
A class that manages a list of Sample objects.
std::vector< Sample * >::const_iterator iterator
the iterator to use
A Sample based on a simple file list.
Definition SampleLocal.h:38
a base class that manages a set of files belonging to a particular data set and the associated meta-d...
Definition Sample.h:54
double entries
Definition listroot.cxx:49
This module provides a lot of global definitions, forward declarations and includes that are used by ...
Definition PrunDriver.h:15
void scanNEvents(SampleHandler &sh)
effects: scan each sample in the sample handler and store the number of entries per file in the meta-...
SampleHandler splitSample(Sample &sample, const Long64_t nevt)
effects: split the given sample into a set of samples, with each sample containing either exactly one...
-diff
static const std::string numEvents
the number of events
Definition MetaFields.h:64
static const std::string treeName_default
the default value of treeName
Definition MetaFields.h:55
static const std::string numEventsPerFile
the number of events in each file
Definition MetaFields.h:67
static const std::string treeName
the name of the tree in the sample
Definition MetaFields.h:52
static const std::string crossSection
the cross section field
Definition MetaFields.h:58
TChain * tree
TFile * file