All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ana::FileReducer Class Reference

Create smaller CAFs. More...

#include <FileReducer.h>

Inheritance diagram for ana::FileReducer:
ana::SpectrumLoaderBase

Public Types

typedef void( ReductionFunc )(caf::StandardRecord *)
 

Public Member Functions

 FileReducer (const std::string &wildcard, const std::string &outfile)
 
 FileReducer (const std::vector< std::string > &fnames, const std::string &outfile)
 
virtual ~FileReducer ()
 
void AddSpillCut (const SpillCut &cut)
 Only copy records to the output file if they pass this cut. More...
 
void AddSliceCut (const SliceCut &cut)
 
void SetEventList (const std::string &fname)
 If called, only events whose run/subrun/event occur in fname will be retained. More...
 
void AddReductionStep (const std::function< ReductionFunc > &f)
 Run the specified reduction function over each event. More...
 
void SetMetadata (const std::string &key, const std::string &val)
 Override any metadata key in the output file. More...
 
virtual void Go () override
 Load all the registered spectra. More...
 

Protected Member Functions

void HandleFile (TFile *fin, TFile *fout, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
 
void HandleNestedTree (TFile *fout, TTree *recTree, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
 
void HandleFlatTree (TFile *fout, TTree *recTree, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
 
void CopyGlobalTree (TFile *fin, TFile *fout)
 
void UpdateMetadata (std::map< std::string, std::string > &meta, const std::set< std::string > &mask, const std::vector< std::string > &fnames) const
 
void Huskify (caf::StandardRecord *sr) const
 Strip all information out of this record and tag it as a husk. More...
 
- Protected Member Functions inherited from ana::SpectrumLoaderBase
 SpectrumLoaderBase (DataSource src=kBeam)
 Component of other constructors. More...
 
 SpectrumLoaderBase (const std::string &wildcard, DataSource src=kBeam)
 Construct from a filename, wildcard, SAM definition, or SAM query. More...
 
 SpectrumLoaderBase (const std::vector< std::string > &fnames, DataSource src=kBeam)
 Construct from an explicit list of files. More...
 
 SpectrumLoaderBase (SpectrumLoaderBase &&)=default
 
SpectrumLoaderBaseoperator= (SpectrumLoaderBase &&)=default
 
 SpectrumLoaderBase (const SpectrumLoaderBase &)=delete
 
SpectrumLoaderBaseoperator= (const SpectrumLoaderBase &)=delete
 
IFileSourceWildcardOrSAMQuery (const std::string &str) const
 Figure out if str is a wildcard or SAM query and return a source. More...
 
virtual void RemoveSpectrum (Spectrum *)
 
virtual void RemoveReweightableSpectrum (ReweightableSpectrum *)
 
int NFiles () const
 Forwards to fFileSource. More...
 
TFile * GetNextFile ()
 
virtual ~SpectrumLoaderBase ()
 
virtual void AddSpectrum (Spectrum &spect, const Var &var, const SpillCut &spillcut, const Cut &cut, const SystShifts &shift, const Var &wei=kUnweighted)
 For use by the Spectrum constructor. More...
 
virtual void AddSpectrum (Spectrum &spect, const MultiVar &var, const SpillCut &spillcut, const Cut &cut, const SystShifts &shift, const Var &wei=kUnweighted)
 For use by the Spectrum constructor. More...
 
virtual void AddSpectrum (Spectrum &spect, const SpillVar &var, const SpillCut &cut, const SpillVar &wei=kSpillUnweighted)
 For use by the Spectrum constructor. More...
 
virtual void AddSpectrum (Spectrum &spect, const SpillMultiVar &var, const SpillCut &cut, const SpillVar &wei=kSpillUnweighted)
 For use by the Spectrum constructor. More...
 
virtual void AddReweightableSpectrum (ReweightableSpectrum &spect, const Var &var, const Cut &cut, const SystShifts &shift, const Var &wei)
 For use by the constructors of ReweightableSpectrum subclasses. More...
 
virtual void AddReweightableSpectrum (ReweightableSpectrum &spect, const Var &var, const SpillCut &spillcut, const SliceCut &slicecut, const SystShifts &shift, const Var &wei)
 For use by the constructors of ReweightableSpectrum subclasses. More...
 
virtual bool Gone () const
 Indicate whether or not Go has been called. More...
 

Protected Attributes

std::string fOutfile
 
SpillCutfSpillCut
 
SliceCutfSliceCut
 
std::set< std::tuple< int, int,
int > > 
fEventList
 
std::vector< std::function
< ReductionFunc > > 
fReductionFuncs
 
std::map< std::string,
std::string > 
fMetaMap
 
- Protected Attributes inherited from ana::SpectrumLoaderBase
std::string fWildcard
 
std::unique_ptr< IFileSourcefFileSource
 
DataSource fSource
 
bool fGone
 Has Go() been called? Can't add more histograms after that. More...
 
double fPOT
 
double fPOTFromHist
 Accumulated by calls to GetNextFile. More...
 
int fNReadouts
 
IDMap< SpillCut, IDMap
< SystShifts, IDMap< Cut,
IDMap< Var, IDMap
< VarOrMultiVar, SpectList > > > > > 
fHistDefs
 All the spectra that need to be filled. More...
 
IDMap< SpillCut, IDMap
< SpillVar, IDMap
< SpillVarOrMultiVar,
SpectList > > > 
fSpillHistDefs
 [spillcut][spillwei][spillvar] More...
 

Additional Inherited Members

- Protected Types inherited from ana::SpectrumLoaderBase
typedef _VarOrMultiVar
< caf::SRSliceProxy
VarOrMultiVar
 
typedef _VarOrMultiVar
< caf::SRSpillProxy
SpillVarOrMultiVar
 

Detailed Description

Create smaller CAFs.

This class produces new CAFs, removing entries that fail a cut. It also allows the event record to be edited in custom ways.

Definition at line 19 of file FileReducer.h.

Member Typedef Documentation

typedef void( ana::FileReducer::ReductionFunc)(caf::StandardRecord *)

Definition at line 36 of file FileReducer.h.

Constructor & Destructor Documentation

ana::FileReducer::FileReducer ( const std::string &  wildcard,
const std::string &  outfile 
)

Definition at line 32 of file FileReducer.cxx.

34  : SpectrumLoaderBase(wildcard),
35  fOutfile(outfile),
36  fSpillCut(nullptr), fSliceCut(nullptr)
37  {
38  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
std::string fOutfile
Definition: FileReducer.h:74
SpectrumLoaderBase(DataSource src=kBeam)
Component of other constructors.
SpillCut * fSpillCut
Definition: FileReducer.h:75
ana::FileReducer::FileReducer ( const std::vector< std::string > &  fnames,
const std::string &  outfile 
)

Definition at line 41 of file FileReducer.cxx.

43  : SpectrumLoaderBase(fnames),
44  fOutfile(outfile),
45  fSpillCut(nullptr),
46  fSliceCut(nullptr)
47  {
48  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
std::string fOutfile
Definition: FileReducer.h:74
SpectrumLoaderBase(DataSource src=kBeam)
Component of other constructors.
SpillCut * fSpillCut
Definition: FileReducer.h:75
ana::FileReducer::~FileReducer ( )
virtual

Definition at line 51 of file FileReducer.cxx.

52  {
53  delete fSpillCut;
54  delete fSliceCut;
55  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
SpillCut * fSpillCut
Definition: FileReducer.h:75

Member Function Documentation

void ana::FileReducer::AddReductionStep ( const std::function< ReductionFunc > &  f)
inline

Run the specified reduction function over each event.

Definition at line 41 of file FileReducer.h.

41 {fReductionFuncs.push_back(f);}
std::vector< std::function< ReductionFunc > > fReductionFuncs
Definition: FileReducer.h:80
void ana::FileReducer::AddSliceCut ( const SliceCut cut)

Definition at line 69 of file FileReducer.cxx.

70  {
71  if(fSliceCut){
72  *fSliceCut = *fSliceCut && cut;
73  }
74  else{
75  fSliceCut = new SliceCut(cut);
76  }
77  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
_Cut< caf::SRSliceProxy > SliceCut
Representation of a cut (selection) to be applied to a caf::StandardRecord object.
Definition: Cut.h:94
void ana::FileReducer::AddSpillCut ( const SpillCut cut)

Only copy records to the output file if they pass this cut.

Definition at line 58 of file FileReducer.cxx.

59  {
60  if(fSpillCut){
61  *fSpillCut = *fSpillCut && cut;
62  }
63  else{
64  fSpillCut = new SpillCut(cut);
65  }
66  }
_Cut< caf::SRSpillProxy > SpillCut
Equivalent of Cut acting on caf::SRSpill. For use in spill-by-spill data quality cuts.
Definition: Cut.h:99
SpillCut * fSpillCut
Definition: FileReducer.h:75
void ana::FileReducer::CopyGlobalTree ( TFile *  fin,
TFile *  fout 
)
protected

Definition at line 322 of file FileReducer.cxx.

323  {
324  TTree* globalIn = (TTree*)fin->Get("globalTree");
325  if(!globalIn) return;
326 
327  // Copy globalTree verbatim from input to output
328  caf::SRGlobal global;
329  caf::SRGlobal* pglobal = &global;
330  globalIn->SetBranchAddress("global", &pglobal);
331  fout->cd();
332  TTree globalOut("globalTree", "globalTree");
333  globalOut.Branch("global", "caf::SRGlobal", &global);
334  assert(globalIn->GetEntries() == 1);
335  // TODO check that the globalTree is the same in every file
336  globalIn->GetEntry(0);
337  globalOut.Fill();
338  globalOut.Write();
339  }
void ana::FileReducer::Go ( )
overridevirtual

Load all the registered spectra.

Implements ana::SpectrumLoaderBase.

Definition at line 95 of file FileReducer.cxx.

96  {
97  // FloatingExceptionOnNaN fpnan;
98 
99  // Don't want overflow to happen. Set to 1 petabyte: effectively infinite.
100  TTree::SetMaxTreeSize(1e15);
101 
102  if(fGone){
103  std::cerr << "Error: can only call Go() once on a FileReducer" << std::endl;
104  abort();
105  }
106  fGone = true;
107 
108  const int Nfiles = NFiles();
109 
110  Progress prog(TString::Format("Filling from %d files matching '%s'", Nfiles, fWildcard.c_str()).Data());
111 
112  TFile fout(fOutfile.c_str(), "RECREATE");
113  TTree* trOut = 0; // created in first call to HandleFile()
114 
115  TH1* hPOTOut = new TH1F("TotalPOT", "", 1, 0, 1);
116  TH1* hEventsOut = new TH1F("TotalEvents", "", 1, 0, 1);
117 
118  std::vector<std::string> fnames;
119 
120  std::map<std::string, std::string> meta;
121  std::set<std::string> meta_mask;
122 
123  long nRecSeen = 0;
124  long nRecPassed = 0;
125 
126  int fileIdx = -1;
127  while(TFile* fin = GetNextFile()){
128  ++fileIdx;
129 
130  assert(!fin->IsZombie());
131 
132  fnames.push_back(fin->GetName());
133 
134  TH1* hPOT = (TH1*)fin->Get("TotalPOT");
135  assert(hPOT);
136  hPOTOut->Add(hPOT);
137 
138  TDirectory* meta_dir = (TDirectory*)fin->Get("metadata");
139  assert(meta_dir);
140  CombineMetadata(meta, GetCAFMetadata(meta_dir), meta_mask);
141 
142  TH1* hEvents = (TH1*)fin->Get("TotalEvents");
143  assert(hEvents);
144  hEventsOut->Add(hEvents);
145 
146  HandleFile(fin, &fout, trOut,
147  Nfiles == 1 ? &prog : 0,
148  nRecSeen, nRecPassed);
149 
150  if(fileIdx == 0) CopyGlobalTree(fin, &fout);
151 
152  prog.SetProgress((fileIdx+1.)/Nfiles);
153  } // end while GetNextFile
154 
155  fout.cd();
156  trOut->Write();
157  hPOTOut->Write();
158  hEventsOut->Write();
159 
160  UpdateMetadata(meta, meta_mask, fnames);
161  WriteCAFMetadata(fout.mkdir("metadata"), meta);
162 
163  fout.Close();
164 
165  prog.Done();
166 
167  std::cout << "Passed " << nRecPassed << " / " << nRecSeen << " records";
168  std::cout << std::endl;
169  }
std::string fOutfile
Definition: FileReducer.h:74
BEGIN_PROLOG could also be cerr
void CopyGlobalTree(TFile *fin, TFile *fout)
bool fGone
Has Go() been called? Can&#39;t add more histograms after that.
void CombineMetadata(std::map< std::string, std::string > &base, const std::map< std::string, std::string > &add, std::set< std::string > &mask)
void HandleFile(TFile *fin, TFile *fout, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
std::map< std::string, std::string > GetCAFMetadata(TDirectory *dir)
void WriteCAFMetadata(TDirectory *dir, const std::map< std::string, std::string > &meta)
void UpdateMetadata(std::map< std::string, std::string > &meta, const std::set< std::string > &mask, const std::vector< std::string > &fnames) const
int NFiles() const
Forwards to fFileSource.
prog
Definition: just_cmake.sh:3
BEGIN_PROLOG could also be cout
void ana::FileReducer::HandleFile ( TFile *  fin,
TFile *  fout,
TTree *&  trOut,
Progress prog,
long &  nRecSeen,
long &  nRecPassed 
)
protected

Definition at line 172 of file FileReducer.cxx.

175  {
176  TTree* recTree = (TTree*)fin->Get("recTree");
177  assert(recTree);
178 
179  const caf::CAFType type = caf::GetCAFType(recTree);
180 
181  if(trOut){
182  const caf::CAFType outtype = caf::GetCAFType(trOut);
183  if(type != outtype){
184  std::cerr << "FileReducer: Error: dataset contains mixed CAF types (flat vs nested)" << std::endl;
185  abort();
186  }
187  }
188 
189  if(type == caf::kNested) HandleNestedTree(fout, recTree, trOut,
190  prog,
191  nRecSeen, nRecPassed);
192  else if(type == caf::kFlat) HandleFlatTree(fout, recTree, trOut,
193  prog,
194  nRecSeen, nRecPassed);
195 
196  else{
197  std::cerr << "FileReducer: Error: Unrecognized file type: "
198  << fin->GetName() << std::endl;
199  abort();
200  }
201  }
BEGIN_PROLOG could also be cerr
void HandleNestedTree(TFile *fout, TTree *recTree, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
void HandleFlatTree(TFile *fout, TTree *recTree, TTree *&trOut, Progress *prog, long &nRecSeen, long &nRecPassed)
prog
Definition: just_cmake.sh:3
void ana::FileReducer::HandleFlatTree ( TFile *  fout,
TTree *  recTree,
TTree *&  trOut,
Progress prog,
long &  nRecSeen,
long &  nRecPassed 
)
protected

Definition at line 282 of file FileReducer.cxx.

285  {
286  if(!fEventList.empty()){
287  std::cerr << "FileReducer: Event list not supported for FlatCAFs (yet)" << std::endl;
288  abort();
289  }
290 
291  if(fSpillCut){
292  std::cerr << "FileReducer: Spill cuts not supported for FlatCAFs" << std::endl;
293  abort();
294  }
295 
296  if(fSliceCut){
297  std::cerr << "FileReducer: Slice cuts not supported for FlatCAFs" << std::endl;
298  abort();
299  }
300 
301  if(!fReductionFuncs.empty()){
302  std::cerr << "FileReducer: Reduction functions not supported for FlatCAFs" << std::endl;
303  abort();
304  }
305 
306  // Do the actual copy
307  if(!trOut){
308  fout->cd();
309  trOut = recTree->CloneTree();
310  }
311  else{
312  fout->cd();
313  recTree->CopyAddresses(trOut);
314  trOut->CopyEntries(recTree);
315  }
316 
317  nRecSeen += recTree->GetEntries();
318  nRecPassed += recTree->GetEntries();
319  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
BEGIN_PROLOG could also be cerr
std::vector< std::function< ReductionFunc > > fReductionFuncs
Definition: FileReducer.h:80
SpillCut * fSpillCut
Definition: FileReducer.h:75
std::set< std::tuple< int, int, int > > fEventList
Definition: FileReducer.h:78
void ana::FileReducer::HandleNestedTree ( TFile *  fout,
TTree *  recTree,
TTree *&  trOut,
Progress prog,
long &  nRecSeen,
long &  nRecPassed 
)
protected

Do we need to include the event? Either based on the selection...

Definition at line 204 of file FileReducer.cxx.

208  {
209  if(!trOut){
210  fout->cd();
211  trOut = new TTree("recTree", "recTree");
212  // FloatingExceptionOnNaN fpnan(false);
213  caf::StandardRecord dummy;
214  trOut->Branch("rec", &dummy);
215  }
216 
217  // Use this one for assessing cuts etc
218  caf::Proxy<caf::StandardRecord> srProxy(recTree, "rec");
219 
220  // And if they pass load into this one for writing out
221  caf::StandardRecord* sr = 0;
222  recTree->SetBranchAddress("rec", &sr);
223 
224  caf::StandardRecord* oldsr = 0;
225 
226  const int Nentries = recTree->GetEntries();
227  for(int n = 0; n < Nentries; ++n){
228  ++nRecSeen;
229  recTree->LoadTree(n);
230 
231  // Apply EventList cut if it's been enabled
232  if(!fEventList.empty() &&
233  !fEventList.count(std::make_tuple(srProxy.hdr.run,
234  srProxy.hdr.subrun,
235  srProxy.hdr.evt))) continue;
236 
237  /// Do we need to include the event? Either based on the selection...
238  const bool passesSpillCut = !fSpillCut || (*fSpillCut)(&srProxy);
239  // Or if it carries POT or spill counting information
240  if(passesSpillCut || srProxy.hdr.pot > 0 || srProxy.hdr.ngenevt > 0){
241  recTree->GetEntry(n);
242 
243  if(sr != oldsr){
244  trOut->SetBranchAddress("rec", &sr);
245  oldsr = sr;
246  }
247 
248  if(!passesSpillCut){
249  // We're only keeping this one for the exposure info
250  Huskify(sr);
251  trOut->Fill();
252  continue;
253  }
254  // Otherwise, need to proceed to keep/reject individual slices
255 
256  std::vector<int> tocut;
257  for(unsigned int i = 0; i < srProxy.slc.size(); ++i){
258  if(fSliceCut && !(*fSliceCut)(&srProxy.slc[i])) tocut.push_back(i);
259  }
260 
261  // Remove slices in reverse order so that the indices remain valid
262  for(auto it = tocut.rbegin(); it != tocut.rend(); ++it){
263  sr->slc.erase(sr->slc.begin() + *it);
264  }
265 
266  // Apply any additional reduction steps
267  for(const auto& f: fReductionFuncs) f(sr);
268  // This is kind of problematic since the proxy and actual record
269  // could be out of sync. Let's just disable this option for now.
270  // for(const auto & f: fReductionFuncsWithProxy) f(sr, &srProxy);
271 
272  ++nRecPassed;
273  trOut->Fill();
274  }
275 
276  // prog won't be passed if the caller doesn't want per-event updates
277  if(n%100 == 0 && prog) prog->SetProgress(double(n)/Nentries);
278  } // end for n
279  }
SliceCut * fSliceCut
Definition: FileReducer.h:76
std::vector< SRSlice > slc
Slice branch.
std::vector< std::function< ReductionFunc > > fReductionFuncs
Definition: FileReducer.h:80
SpillCut * fSpillCut
Definition: FileReducer.h:75
The StandardRecord is the primary top-level object in the Common Analysis File trees.
std::set< std::tuple< int, int, int > > fEventList
Definition: FileReducer.h:78
void Huskify(caf::StandardRecord *sr) const
Strip all information out of this record and tag it as a husk.
prog
Definition: just_cmake.sh:3
void ana::FileReducer::Huskify ( caf::StandardRecord sr) const
protected

Strip all information out of this record and tag it as a husk.

Definition at line 384 of file FileReducer.cxx.

385  {
386  sr->hdr.husk = true;
387 
388  // It's not actually necessary for correctness to clear all these fields
389  // out, so not a disaster if we happen to miss one that gets added in
390  // future, for example. But it does help to minimize the size of the
391  // resulting CAF.
392 
393  sr->slc.clear();
394  sr->nslc = 0;
395  sr->fake_reco.clear();
396  sr->nfake_reco = 0;
397  sr->true_particles.clear();
398  sr->ntrue_particles = 0;
399  sr->crt_hits.clear();
400  sr->ncrt_hits = 0;
401  sr->crt_tracks.clear();
402  sr->ncrt_tracks = 0;
403 
404  sr->reco.trk.clear();
405  sr->reco.ntrk = 0;
406  sr->reco.shw.clear();
407  sr->reco.nshw = 0;
408  sr->reco.stub.clear();
409  sr->reco.nstub = 0;
410 
411  sr->mc.nu.clear();
412  sr->mc.nnu = 0;
413  sr->mc.prtl.clear();
414  sr->mc.nprtl = 0;
415  }
int ncrt_tracks
Number of CRT tracks in event.
SRHeader hdr
Header branch: run, subrun, etc.
std::vector< SRMeVPrtl > prtl
If present – information on decay of MeV &quot;Portal&quot; particle.
Definition: SRTruthBranch.h:24
int nfake_reco
Number of Fake-Reco&#39;s in list.
int nslc
Number of slices in list.
std::vector< SRShower > shw
Vector of trac showers.
std::vector< SRSlice > slc
Slice branch.
size_t nshw
Number of trac showers.
size_t nprtl
Number of portals.
Definition: SRTruthBranch.h:25
int ncrt_hits
Number of CRT hits in event.
std::vector< SRTrueInteraction > nu
Vector of true nu or cosmic.
Definition: SRTruthBranch.h:21
size_t nstub
Number of stubs.
std::vector< SRCRTHit > crt_hits
CRT hits in event.
size_t ntrk
Number of panora tracks.
std::vector< SRTrueParticle > true_particles
True particles in spill.
std::vector< SRTrack > trk
Vector of pandora tracks.
std::vector< SRCRTTrack > crt_tracks
CRT tracks in event.
int ntrue_particles
Number of true particles in list.
std::vector< SRFakeReco > fake_reco
List of fake-reco slices.
size_t nnu
Number of true nu or cosmic.
Definition: SRTruthBranch.h:22
SRSliceRecoBranch reco
Slice reco branch: tracks, showers, etc.
SRTruthBranch mc
Truth branch for all interactions.
std::vector< SRStub > stub
Vector of stubs.
void ana::FileReducer::SetEventList ( const std::string &  fname)

If called, only events whose run/subrun/event occur in fname will be retained.

Definition at line 80 of file FileReducer.cxx.

81  {
82  FILE* f = fopen(fname.c_str(), "r");
83  assert(f);
84 
85  while(!feof(f)){
86  int run, subrun, event;
87  fscanf(f, "%d %d %d", &run, &subrun, &event);
88  fEventList.emplace(run, subrun, event);
89  }
90 
91  fclose(f);
92  }
string fname
Definition: demo.py:5
std::set< std::tuple< int, int, int > > fEventList
Definition: FileReducer.h:78
void ana::FileReducer::SetMetadata ( const std::string &  key,
const std::string &  val 
)
inline

Override any metadata key in the output file.

Definition at line 45 of file FileReducer.h.

46  {
47  // Let's assume the user is just setting string keys
48  fMetaMap[key] = "\""+val+"\"";
49  }
std::map< std::string, std::string > fMetaMap
Definition: FileReducer.h:83
void ana::FileReducer::UpdateMetadata ( std::map< std::string, std::string > &  meta,
const std::set< std::string > &  mask,
const std::vector< std::string > &  fnames 
) const
protected

Definition at line 342 of file FileReducer.cxx.

345  {
346  for(const std::string& m: mask){
347  std::cerr << "Warning: metadata parameter '" << m << "' differs between input files and has been dropped from the output." << std::endl;
348  meta.erase(m);
349  }
350 
351  /*
352  // change caf -> decaf in the metadata field, if parents have data_tier
353  if(meta.find("data_tier") != meta.end()){
354  std::string decaf_tier = meta["data_tier"];
355  assert(decaf_tier.size() >= 3);
356  // For indexing: note that all of this is surrounded by quotes
357  assert(decaf_tier.substr(decaf_tier.size()-4,3) == "caf");
358  // don't make 'decaf' into 'dedecaf', however
359  if (decaf_tier.size() < 6 || decaf_tier.substr(decaf_tier.size()-6,5) != "decaf")
360  decaf_tier.replace(decaf_tier.size()-4,3,"decaf");
361  meta["data_tier"] = decaf_tier;
362  }
363  */
364 
365  std::string parents = "[";
366  for(const std::string& f: fnames){
367  parents += "{\"file_name\": \""+std::string(basename((char *)f.c_str()))+"\"}, ";
368  }
369  if(!fnames.empty()) parents.resize(parents.size()-2); // drop trailing ", "
370  parents += "]";
371 
372  meta["parents"] = parents;
373 
374  // if there's one more than one parent, this is a concat
375  if(fnames.size() > 1 && meta.find("data_tier") != meta.end()){
376  meta["data_tier"] = "\"concat_caf\"";
377  }
378 
379  // Allow user to override any metadata
380  for(auto it: fMetaMap) meta[it.first] = it.second;
381  }
BEGIN_PROLOG could also be cerr
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
std::map< std::string, std::string > fMetaMap
Definition: FileReducer.h:83

Member Data Documentation

std::set<std::tuple<int, int, int> > ana::FileReducer::fEventList
protected

Definition at line 78 of file FileReducer.h.

std::map<std::string, std::string> ana::FileReducer::fMetaMap
protected

Definition at line 83 of file FileReducer.h.

std::string ana::FileReducer::fOutfile
protected

Definition at line 74 of file FileReducer.h.

std::vector<std::function<ReductionFunc> > ana::FileReducer::fReductionFuncs
protected

Definition at line 80 of file FileReducer.h.

SliceCut* ana::FileReducer::fSliceCut
protected

Definition at line 76 of file FileReducer.h.

SpillCut* ana::FileReducer::fSpillCut
protected

Definition at line 75 of file FileReducer.h.


The documentation for this class was generated from the following files: