ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/src/AlignmentIORootBase.cc
Revision: 1.4
Committed: Thu Nov 30 09:52:49 2006 UTC (18 years, 5 months ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_1_3_6, CMSSW_1_3_1_HLT6, CMSSW_1_3_5, CMSSW_1_3_1_HLT5, CMSSW_1_3_1_HLT4, CMSSW_1_3_1_HLT3, CMSSW_1_3_4, CMSSW_1_3_3, CMSSW_1_3_2, CMSSW_1_3_1, CMSSW_1_4_0_pre2, CMSSW_1_3_0, CMSSW_1_3_0_pre7, CMSSW_1_4_0_pre1, CMSSW_1_3_0_pre6, cklae_20070314, V01-02-00, V01-01-04, CMSSW_1_3_0_pre5, CMSSW_1_3_0_pre4, V01-01-03, CMSSW_1_3_0_pre3, V01-01-02, CMSSW_1_3_0_SLC4_pre2, CMSSW_1_3_0_pre2, V01-01-01, V01-01-00, CMSSW_1_3_0_SLC4_pre1, CMSSW_1_3_0_pre1, V01-00-00, V00-08-01, V00-08-00
Changes since 1.3: +1 -1 lines
Log Message:
change tree names from <name>:<n> to <name>_<n> to allow easy drawing with
tree-friends

File Contents

# Content
1 // this class's header
2 #include "FWCore/MessageLogger/interface/MessageLogger.h"
3
4 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentIORootBase.h"
5
6 // ----------------------------------------------------------------------------
7 // open file/trees for write
8
9 int AlignmentIORootBase::openRoot(const char* filename, int iteration, bool write)
10 {
11 bWrite=write;
12 int iter;
13
14 edm::LogInfo("AlignmentIORootBase") << "File: " << filename ;
15
16 if (bWrite) { // writing
17
18 int iterfile = testFile(filename,treename);
19 if (iterfile == -1) {
20 iter=iteration;
21 edm::LogInfo("AlignmentIORootBase") << "Write to new file; first iteration: " << iter ;
22 IORoot = new TFile(filename,"recreate");
23 } else {
24 if (iteration == -1) {
25 iter=iterfile+1;
26 edm::LogInfo("AlignmentIORootBase")
27 << "Write to existing file; highest iteration: " << iter;
28 }
29 else {
30 if (iteration<=iterfile) {
31 edm::LogError("AlignmentIORootBase")
32 << "Iteration " << iteration
33 <<" invalid or already exists for tree " << treename;
34 return -1;
35 }
36 iter = iteration;
37 edm::LogInfo("AlignmentIORootBase") << "Write to new iteration: " << iter;
38 }
39 IORoot = new TFile(filename,"update");
40
41 }
42
43 IORoot->cd();
44
45 // create tree
46 edm::LogInfo("AlignmentIORootBase") << "Tree: " << treeName(iter,treename);
47 tree = new TTree(treeName(iter,treename),treetxt);
48 createBranches();
49
50 } else { // reading
51
52 int iterfile = testFile(filename,treename);
53 if ( iterfile == -1 ) {
54 edm::LogError("AlignmentIORootBase") << "File does not exist!";
55 return -1;
56 } else if ( iterfile == -2 ) {
57 edm::LogError("AlignmentIORootBase") << "Tree " << treename
58 << " does not exist in file " << filename;
59 return -1;
60 } else {
61 if (iteration == -1) {
62 iter=iterfile;
63 edm::LogInfo("AlignmentIORootBase") << "Read from highest iteration: " << iter;
64 } else {
65 if (iteration>iterfile) {
66 edm::LogError("AlignmentIORootBase")
67 <<"Iteration " << iteration << " does not exist for tree " << treename;
68 return -1;
69 }
70 iter = iteration;
71 edm::LogInfo("AlignmentIORootBase") << "Read from specified iteration: " << iter;
72 }
73 IORoot = new TFile(filename,"update");
74 }
75
76 IORoot->cd();
77 // set trees
78 edm::LogInfo("AlignmentIORootBase") <<" Tree: " <<treeName(iter,treename);
79 tree = (TTree*)IORoot->Get(treeName(iter,treename));
80 if (tree==NULL) {
81 edm::LogError("AlignmentIORootBase") <<"Tree does not exist in file!";
82 return -1;
83 }
84 setBranchAddresses();
85 }
86
87 return 0;
88 }
89
90 // ----------------------------------------------------------------------------
91 // write tree and close file
92
93 int AlignmentIORootBase::closeRoot(void)
94 {
95 if (bWrite) { //writing
96 IORoot->cd();
97 tree->Write();
98 tree->Delete();
99 IORoot->Close();
100 }
101 else { // reading
102 IORoot->cd();
103 tree->Delete();
104 IORoot->Close();
105 }
106
107 return 0;
108 }
109
110 // ----------------------------------------------------------------------------
111 // returns highest existing iteration in file
112 // if file does not exist: return -1
113
114 int AlignmentIORootBase::testFile(const char* filename, const TString &tname)
115 {
116 FILE* testFILE;
117 testFILE = fopen(filename,"r");
118 if (testFILE == NULL) {
119 return -1;
120 }
121 else {
122 fclose(testFILE);
123 int ihighest=-2;
124 TFile IORoot(filename,"update");
125 for (int iter=0; iter<itermax; iter++) {
126 if ((0 != (TTree*)IORoot.Get(treeName(iter,tname)))
127 && (iter>ihighest)) ihighest=iter;
128 }
129 IORoot.Close();
130 return ihighest;
131 }
132 }
133
134 // ----------------------------------------------------------------------------
135 // create tree name from stub+iteration
136
137 TString AlignmentIORootBase::treeName(int iter, const TString &tname)
138 {
139 return TString(tname + Form("_%i",iter));
140 }