ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeService/src/TreeService.cc
Revision: 1.6
Committed: Fri Jun 20 17:52:24 2008 UTC (16 years, 10 months ago) by loizides
Content type: text/plain
Branch: MAIN
Changes since 1.5: +38 -9 lines
Log Message:
Check for existence of parameters. This allows one to specify the TreeService in almost a single line since now meaningful default values can be set.

File Contents

# User Rev Content
1 loizides 1.6 // $Id: TreeService.cc,v 1.5 2008/06/18 13:23:23 paus Exp $
2 loizides 1.1
3     #include "MitProd/TreeService/interface/TreeService.h"
4     #include "DataFormats/Provenance/interface/ModuleDescription.h"
5     #include "FWCore/ParameterSet/interface/ParameterSet.h"
6     #include "FWCore/ServiceRegistry/interface/ActivityRegistry.h"
7     #include "FWCore/ServiceRegistry/interface/Service.h"
8     #include "FWCore/MessageLogger/interface/JobReport.h"
9     #include "MitAna/DataUtil/interface/TreeWriter.h"
10 loizides 1.6 #include "MitAna/DataTree/interface/Names.h"
11 loizides 1.1
12     using namespace edm;
13     using namespace std;
14     using namespace mithep;
15    
16 paus 1.2 //--------------------------------------------------------------------------------------------------
17 loizides 1.3 TreeService::TreeService(const ParameterSet &cfg, ActivityRegistry &r) :
18 loizides 1.6 tws_ (0)
19 loizides 1.3 {
20 loizides 1.6 if (cfg.exists("treeNames"))
21     treeNames_=cfg.getUntrackedParameter<vector<string> >("treeNames");
22     else
23     treeNames_.push_back(Names::gkEvtTreeName);
24    
25     if (cfg.exists("fileNames"))
26     fileNames_=cfg.getUntrackedParameter<vector<string> >("fileNames");
27     else
28     fileNames_.push_back("mit-test");
29    
30     if (cfg.exists("pathNames"))
31     pathNames_=cfg.getUntrackedParameter<vector<string> >("pathNames");
32     else
33     pathNames_.push_back(".");
34    
35     if (cfg.exists("maxSizes"))
36     maxSizes_=cfg.getUntrackedParameter<vector<unsigned> >("maxSizes");
37     else
38     maxSizes_.push_back(1024);
39    
40     if (cfg.exists("compLevels"))
41     compLevels_=cfg.getUntrackedParameter<vector<unsigned> >("compLevels");
42     else
43     compLevels_.push_back(9);
44    
45     if (cfg.exists("splitLevels"))
46     splitLevels_=cfg.getUntrackedParameter<vector<unsigned> >("splitLevels");
47     else
48     splitLevels_.push_back(99);
49    
50     if (cfg.exists("brSizes"))
51     brSizes_=cfg.getUntrackedParameter<vector<unsigned> >("brSizes");
52     else
53     brSizes_.push_back(32*1024);
54    
55 loizides 1.3 if (treeNames_.size()!=fileNames_.size()) {
56     throw edm::Exception(edm::errors::Configuration, "TreeService::TreeService()\n")
57     << "Number of main trees should match number of files " << treeNames_.size()
58     << " " << fileNames_.size() << "\n";
59     return;
60 loizides 1.1 }
61    
62     // setup tw array
63     tws_.SetOwner(kTRUE);
64     tws_.Expand(treeNames_.size());
65    
66     // init tree writers
67 paus 1.2 for (unsigned int i=0; i<treeNames_.size(); ++i) {
68 loizides 1.1
69     TreeWriter *t = new TreeWriter(treeNames_.at(i).c_str(),1);
70    
71     t->SetPrefix(fileNames_.at(i).c_str());
72 paus 1.2
73     if (i<pathNames_.size())
74     t->SetBaseURL(pathNames_.at(i).c_str());
75     else if (pathNames_.size()>0)
76     t->SetBaseURL(pathNames_.at(0).c_str());
77    
78     if (i<maxSizes_.size())
79 loizides 1.4 t->SetMaxSize((Long64_t)maxSizes_.at(i)*1024*1024);
80 paus 1.2 else if (maxSizes_.size()>0)
81 loizides 1.4 t->SetMaxSize((Long64_t)maxSizes_.at(0)*1024*1024);
82 paus 1.2
83     if (i<compLevels_.size())
84     t->SetCompressLevel(compLevels_.at(i));
85     else if (compLevels_.size()>0)
86     t->SetCompressLevel(compLevels_.at(0));
87    
88     if (i<splitLevels_.size())
89     t->SetDefaultSL(splitLevels_.at(i));
90     else if (splitLevels_.size()>0)
91     t->SetDefaultSL(splitLevels_.at(0));
92    
93     if (i<brSizes_.size())
94     t->SetDefaultBrSize(brSizes_.at(i));
95     else if (brSizes_.size()>0)
96     t->SetDefaultBrSize(brSizes_.at(0));
97 loizides 1.1
98     //t->Print();
99     tws_.Add(t);
100     }
101    
102     // set watchers
103 paus 1.2 r.watchPreProcessEvent (this,&TreeService::preEventProcessing);
104 loizides 1.1 r.watchPostProcessEvent(this,&TreeService::postEventProcessing);
105 paus 1.2 r.watchPostBeginJob (this,&TreeService::postBeginJob);
106     r.watchPostEndJob (this,&TreeService::postEndJob);
107 loizides 1.1 }
108    
109 paus 1.2 //--------------------------------------------------------------------------------------------------
110 loizides 1.1 TreeService::~TreeService()
111     {
112     }
113    
114 paus 1.2 //--------------------------------------------------------------------------------------------------
115 loizides 1.1 TreeWriter* TreeService::get(const char *name)
116     {
117 paus 1.2 if (tws_.GetEntries()<=0)
118     return 0;
119 loizides 1.1
120 paus 1.2 if (name==0)
121 loizides 1.1 return dynamic_cast<TreeWriter*>(tws_.At(0));
122    
123     TObject *ret = tws_.FindObject(name);
124 paus 1.2
125 loizides 1.1 return dynamic_cast<TreeWriter*>(ret);
126     }
127    
128 paus 1.2 //--------------------------------------------------------------------------------------------------
129 loizides 1.1 void TreeService::postBeginJob()
130     {
131     // nothing to be done for now
132     }
133    
134 paus 1.2 //--------------------------------------------------------------------------------------------------
135 loizides 1.1 void TreeService::postEndJob()
136     {
137     tws_.Clear();
138     }
139    
140 paus 1.2 //--------------------------------------------------------------------------------------------------
141 loizides 1.1 void TreeService::postEventProcessing(const Event&, const EventSetup&)
142     {
143 paus 1.2 for (int i=0; i<tws_.GetEntries(); ++i) {
144 loizides 1.1 TreeWriter *tw=dynamic_cast<TreeWriter*>(tws_.At(i));
145 paus 1.2 if (tw)
146 loizides 1.1 tw->EndEvent();
147     }
148     }
149    
150 paus 1.2 //--------------------------------------------------------------------------------------------------
151 loizides 1.1 void TreeService::preEventProcessing(const EventID&, const Timestamp&)
152     {
153 paus 1.2 for (int i=0; i<tws_.GetEntries(); ++i) {
154 loizides 1.1 TreeWriter *tw=dynamic_cast<TreeWriter*>(tws_.At(i));
155 paus 1.2 if (tw)
156 loizides 1.1 tw->BeginEvent();
157     }
158     }