BOSS 6.6.4.p03
BESIII Offline Software System
Loading...
Searching...
No Matches
RootInterface.cxx
Go to the documentation of this file.
1#include "TFile.h"
2#include "TTree.h"
3#include "TBranch.h"
4#include "TChain.h"
5#include "TClonesArray.h"
6
11
12#include "GaudiKernel/SmartIF.h"
13#include "GaudiKernel/Bootstrap.h"
14#include "GaudiKernel/IAppMgrUI.h"
15#include "GaudiKernel/IProperty.h"
16#include "GaudiKernel/ISvcLocator.h"
17
18#include <algorithm> // for find
19#include <iostream>
20#include <fstream>
21#include <string>
22RootInterface* RootInterface::m_rootInterface = 0;
23
24
26 if (m_rootInterface) return m_rootInterface;
27 m_rootInterface=new RootInterface(log);
28 return m_rootInterface;
29}
30
31RootInterface::RootInterface(MsgStream str) : log(str)
32{
33
34 m_branches= new TClonesArray("TBranch",1);
35 m_branchesRead= new TClonesArray("TBranch",1);
36 m_entries=-1;
37 m_EOF=false;
38 m_ENDFILE = false ;
39 m_fileNum = 0; //-1
40 }
41
42
44{
45}
46
48{
49 IInterface* iface = Gaudi::createApplicationMgr();
50 SmartIF<IProperty> propMgr (iface );
51 std::string path;
52 propMgr->getProperty( "JobOptionsPath", path);
53 log << MSG::INFO << "JobOptions file for current job: " <<path << endreq;
54 ifstream fin(path.c_str());
55 string jobOptions;
56 string tempString;
57 while(getline(fin,tempString))
58 {
59 if( tempString.size()>0 && tempString.find("//")>tempString.size() )
60 {
61 jobOptions += tempString;
62 jobOptions += "\n";
63 }
64 }
65 log << MSG::INFO << "JobOptions: " << endreq
66 << jobOptions << endreq;
67 return jobOptions;
68}
69
71{
72 ISvcLocator* svcLocator = Gaudi::svcLocator();
73 IDataInfoSvc *tmpInfoSvc;
74 DataInfoSvc* jobInfoSvc;
75 string decayOptions;
76 StatusCode status = svcLocator->service("DataInfoSvc",tmpInfoSvc);
77 if (status.isSuccess()) {
78 log << MSG::INFO << "get the DataInfoSvc" << endreq;
79 jobInfoSvc=dynamic_cast<DataInfoSvc *>(tmpInfoSvc);
80 decayOptions = jobInfoSvc->getDecayOptions();
81 log << MSG::INFO << "get decay options" << endreq
82 << decayOptions << endreq;
83 }else {
84 log << MSG::WARNING << "could not get the DataInfoSvc. Ignore it." << endreq;
85 }
86 return decayOptions;
87}
88
90{
91 ISvcLocator* svcLocator = Gaudi::svcLocator();
92 IDataInfoSvc *tmpInfoSvc;
93 DataInfoSvc* jobInfoSvc;
94 std::vector<int> totEvtNo;
95 StatusCode status = svcLocator->service("DataInfoSvc",tmpInfoSvc);
96 if (status.isSuccess()) {
97 log << MSG::INFO << "get the DataInfoSvc" << endreq;
98 jobInfoSvc=dynamic_cast<DataInfoSvc *>(tmpInfoSvc);
99 totEvtNo = jobInfoSvc->getTotEvtNo();
100 log << MSG::INFO << "get total event number for each run" << endreq;
101 }else {
102 log << MSG::WARNING << "could not get the DataInfoSvc. Ignore it." << endreq;
103 }
104 return totEvtNo;
105}
106
108
109 // Get the messaging service, print where you are
110 log << MSG::INFO << "finalize() in RootInterface" << endreq;
111
112 // close file (FIXME for several output files)
113 std::vector<TTree *>::const_iterator trees;
114 for (trees=m_outputTrees.begin();trees<m_outputTrees.end();trees++)
115 if (*trees) {
116 int treenr=(*trees)->GetUniqueID();
117 if (m_outputFiles[treenr]) {
118 if (!m_outputFiles[treenr]->IsOpen()) {
119 log << MSG::ERROR << "Could not open file for writing" << endreq;
120 return StatusCode::FAILURE;
121 } else {
122 log << MSG::DEBUG<<" Closing file "<<treenr<<", tree "<<(*trees)->GetName()<<endreq;
123 TDirectory *saveDir = gDirectory;
124 m_outputFiles[treenr]->cd();
125
126 TJobInfo* jobInfo = new TJobInfo;
127 TTree* m_jobInfoTree = new TTree("JobInfoTree","Job info");
128 m_jobInfoTree->Branch("JobInfo",&jobInfo);
129
130 m_bossVer = getenv("BES_RELEASE");
131 log << MSG::INFO << "fill boss version: "<<m_bossVer << endreq;
132
133 string tmpJobOptions = getJobOptions();
134 m_jobOptions.push_back( tmpJobOptions );
135
136 if(m_decayOptions.size()==0)
137 m_decayOptions = getDecayOptions();
138
139 m_totEvtNo = getTotEvtNo();
140 jobInfo->setBossVer(m_bossVer);
141 jobInfo->setJobOptions(m_jobOptions);
142 jobInfo->setDecayOptions(m_decayOptions);
143 jobInfo->setTotEvtNo(m_totEvtNo);
144 m_jobInfoTree->Fill();
145
146 //? m_mcFile->Write(0, TObject::kOverwrite);
147 int st =1;
148 st = m_outputFiles[treenr]->Write();
149 if(st==0){
150 log << MSG::FATAL<<" can not write the file "<< m_outputFilenames[treenr].c_str()<<endreq;
151 exit(1);
152 }
153
154 m_outputFiles[treenr]->Close();
155 saveDir->cd();
156 }
157 }
158 }
159 if(m_outputTrees.size()>0) m_outputTrees.clear();
160 return StatusCode::SUCCESS;
161}
162
163StatusCode RootInterface::addInput(const std::string& treename,const std::string& file) {
164 log << MSG::DEBUG << "addInput for Tree "<<treename<<endreq;
165 StatusCode sc=StatusCode::SUCCESS;
166 m_fileNames.push_back(file); //input files 2005-11-28
167 m_otherTrees.push_back(NULL);
168 inputFiles.push_back(NULL);
169 unsigned int treenr;
170 sc=getTreeNr(treename,treenr,true);
171 m_inputFilenames[treenr]=file; // the last one file is setted
172 m_inputFiles[treenr]=NULL;
173 m_inputTrees[treenr]=NULL;
174
175 return sc;
176}
177
178StatusCode RootInterface::addOutput(const std::string& treename,const std::string& file,int split, int bufsize, int compression) {
179 static int i =0;
180 i ++;
181 log << MSG::DEBUG << "addOutput for Tree "<<treename<<endreq;
182 StatusCode sc=StatusCode::SUCCESS;
183 unsigned int treenr;
184 sc=getTreeNr(treename,treenr,true);
185 m_outputFilenames[treenr]=file;
186 m_outputFiles[treenr]=NULL;
187 m_outputTrees[treenr]=NULL;
188 m_splitModes[treenr]=split;
189 m_bufSizes[treenr]=bufsize;
190 m_compressionLevels[treenr]=compression;
191
192 return sc;
193}
194
195StatusCode RootInterface::createBranch(const std::string& treename,const std::string& branchname ,const char *classname,void *addr,int & branchnr){
196
197 log << MSG::DEBUG << "CreateBranch, Tree "<<treename<<" branch "<<branchname<<endreq;
198
199 TBranch *branch;
200 unsigned int treenr;
201 StatusCode sc=getTreeNr(treename,treenr);
202 if (!sc.isSuccess()) return sc;
203
204 if ( m_outputFilenames[treenr].empty() ) {
205 log << MSG::DEBUG << "No corresponding output file specified, ignore createBranch: "<<branchname<<endreq;
206 return StatusCode::SUCCESS;
207 }
208
209 if(!m_outputTrees[treenr]) sc=this->createTree(treenr,treename);
210 if (!sc.isSuccess()) return sc;
211 TTree * tree=m_outputTrees[treenr];
212 tree->SetUniqueID(treenr);
213
214 branch = tree->Branch(branchname.c_str(),classname,addr,m_bufSizes[treenr],m_splitModes[treenr]);
215 branch->SetUniqueID(treenr);
216 branchnr=m_branches->GetEntriesFast()+1;
217 m_branches->Expand(branchnr);
218 TClonesArray &a = *m_branches;
219 a[branchnr-1]=branch;
220 tree->SetBasketSize(branchname.c_str(),m_bufSizes[treenr]); //some problem with above method to set buffersize, so we set it here.
221 return StatusCode::SUCCESS;
222}
223
224
225StatusCode RootInterface::createTree(const unsigned int treenr,const std::string treename)
226{
227 // opens file and creates TTree on it
228
229 TDirectory *saveDir = gDirectory;
230
231 // Create the new ROOT file
232 m_outputFiles[treenr] =new TFile(m_outputFilenames[treenr].c_str(), "RECREATE");
233 if(m_outputFiles[treenr]->IsZombie()){
234 std::cout<<"RootInterface ERROR::Can't not open file"<<m_outputFilenames[treenr].c_str()<<std::endl;
235 exit(1);
236 }
237 if (!m_outputFiles[treenr]->IsOpen()) {
238 log << MSG::FATAL << "ROOT file " << m_outputFilenames[treenr]
239 << " could not be opened for writing." << endreq;
240 exit(1);
241 return StatusCode::FAILURE;
242 }
243 log << MSG::INFO << "RootInterface::opened file for output:" << m_outputFilenames[treenr].c_str() <<endreq;
244
245 m_outputFiles[treenr]->cd();
246 m_outputFiles[treenr]->SetCompressionLevel(m_compressionLevels[treenr]);
247 std::string title=treename+" from conversion";
248 m_outputTrees[treenr]= new TTree(treename.c_str(), title.c_str());
249 TTree::SetMaxTreeSize(20000000000LL);
250
251 saveDir->cd();
252
253 return StatusCode::SUCCESS;
254}
255
256TTree * RootInterface::getTree(const std::string treename)
257{
258
259 // get TTree for input
260 log<<MSG::INFO<<"RootInterface:;getTree"<<endreq;
261 unsigned int treenr;
262 getTreeNr(treename,treenr);
263
264 if (m_inputTrees[treenr]) return m_inputTrees[treenr];
265 if (!m_inputFiles[treenr] ) {
266 m_inputFiles[treenr] = TFile::Open(m_fileNames[treenr].c_str(), "READ");
267 if (!m_inputFiles[treenr]->IsOpen()) {
268 log << MSG::ERROR << "ROOT file " << m_inputFiles[treenr]->GetName()
269 << " could not be opened for reading." << endreq;
270 delete m_inputFiles[treenr];
271 m_inputFiles[treenr]=NULL;
272 m_EOF=true;
273 return NULL;
274 }
275 }
276 log << MSG::INFO << "RootInterface::opened file for input:" << m_fileNames[treenr].c_str() <<endreq;
277 m_currentFileName = m_fileNames[treenr].c_str(); //liangyt 2008-11-19
278 TTree *tree= (TTree *)m_inputFiles[treenr]->Get(treename.c_str());
279 if (!tree) {
280 log << MSG::ERROR << "ROOT file " << m_inputFiles[treenr]->GetName()
281 << " does not contain requested TTree: " << treename <<endreq;
282 return NULL;
283 }
284 if (tree->GetEntries()<=0)
285 {
286 log << MSG::ERROR << "ROOT file "<< m_inputFiles[treenr]->GetName()<< " entries <= 0" << endreq;
287 exit(1);
288 }
289
290 m_inputTrees[treenr]=tree;
291 if (m_entries<=0 ) {
292 m_entries=(Int_t)tree->GetEntries();
293 }
294
295 printJobInfo( m_inputFiles[treenr],1);
296
297 return tree;
298}
299
300void RootInterface::printJobInfo(TFile* file, int level)
301{
302 TTree* tree2 = (TTree *)file->Get("JobInfoTree");
303 if(!tree2)
304 {
305 std::cout<<"no JobInfoTree for file "<<file->GetName()<<std::endl;
306 exit(1);
307 }
308 else
309 {
310 log << MSG::INFO << "get JobInfoTree" << endreq;
311 TBranch* branch = tree2->GetBranch("JobInfo");
312 if(!branch)
313 {
314 std::cout<<"ERROR! No branch in JobInfoTree"<<std::endl;
315 exit(1);
316 }
317 else
318 {
319 TJobInfo* jobInfo= new TJobInfo;
320 branch->SetAddress(&jobInfo);
321 branch->GetEntry(0);
322 m_bossVer = jobInfo->getBossVer() ;
323 std::cout<<std::endl
324 << "**************************************************" << std::endl
325 <<"Print JobInfo for data file: "<< file->GetName()<<std::endl
326 << " BOSS version: "<< m_bossVer << std::endl
327 << "**************************************************" << std::endl
328 << std::endl;
329
330 m_decayOptions = jobInfo->getDecayOptions();
331 if(m_decayOptions.size()>0)
332 {
333 std::cout<< std::endl
334 <<"**************************************************" << std::endl
335 <<" Decay Options: "<<std::endl
336 <<m_decayOptions << std::endl
337 <<"**************************************************" << std::endl
338 << std::endl;
339 }
340
341 ISvcLocator* svcLocator = Gaudi::svcLocator();
342 IDataInfoSvc *tmpInfoSvc;
343 DataInfoSvc* jobInfoSvc;
344 StatusCode status = svcLocator->service("DataInfoSvc",tmpInfoSvc);
345 if (status.isSuccess()) {
346 log << MSG::INFO << "get the DataInfoSvc" << endreq;
347 jobInfoSvc=dynamic_cast<DataInfoSvc *>(tmpInfoSvc);
348 }else {
349 log << MSG::WARNING << "could not get the DataInfoSvc." << endreq;
350 }
351
352 m_totEvtNo = jobInfo->getTotEvtNo();
353 jobInfoSvc->setTotEvtNo(m_totEvtNo);
354
355 if(level>0)
356 {
357 std::cout<<std::endl
358 <<"**************************************************" << std::endl
359 <<" JobOptions for this data file: " << std::endl
360 << std::endl;
361
362
363 m_jobOptions = jobInfo->getJobOptions();
364 vector<std::string> vs = m_jobOptions;
365 int nv = vs.size();
366 if(nv>0)
367 {
368 for(int i=0;i<nv;i++)
369 {
370 std::cout<<vs[i]<<std::endl;
371 std::cout<<" end of the jobOptions file " << std::endl;
372 std::cout<<"**************************************************" << std::endl
373 <<std::endl;
374 }
375 }
376 }
377 }
378 }
379}
380
381TTree * RootInterface::getOtherTree(const std::string treename){
382 //get other TTree for input
383 log<<MSG::INFO<<"RootInterface:;getOtherTree"<<endreq;
384 m_ENDFILE = false ;
385 if(m_otherTrees[m_fileNum]) return m_otherTrees[m_fileNum];
386 // TFile* inputFile = new TFile(m_fileNames[m_fileNum].c_str(),"READ");
387 inputFiles[m_fileNum] = TFile::Open(m_fileNames[m_fileNum].c_str(),"READ");
388
389 if(!inputFiles[m_fileNum]->IsOpen()){
390 log<<MSG::ERROR<<"ROOT File" <<inputFiles[m_fileNum]->GetName()<<"Coult not be opened for reading."<<endreq;
391 delete inputFiles[m_fileNum];
392 inputFiles[m_fileNum] = NULL;
393 return NULL; // The Root can not be opened
394 }
395 m_EOF = false;
396 log<<MSG::INFO<<"RootIntrFace:;Opened File for input:"<<m_fileNames[m_fileNum].c_str()<<endreq;
397 m_currentFileName = m_fileNames[m_fileNum].c_str(); //liangyt 2008-11-19
398
399 TTree* tree =(TTree*)inputFiles[m_fileNum]->Get(treename.c_str());//the same tree name;
400 if(!tree){
401 log << MSG::ERROR << "ROOT file " << inputFiles[m_fileNum]->GetName()
402 << " does not contain requested TTree: " << treename <<endreq;
403 return NULL;
404 }
405
406 if(tree->GetEntries()<=0)
407 {
408 log << MSG::ERROR << "ROOT file "<< m_fileNames[m_fileNum].c_str()<< " entries <= 0" << endreq;
409 exit(1);
410 }
411
412 m_otherTrees[m_fileNum] = tree;
413 if (m_entries<=0){
414 m_entries=(Int_t)tree->GetEntries();
415 log<<MSG::INFO<<"m_entries = "<<m_entries<<endreq;
416 }
417
418 printJobInfo(inputFiles[m_fileNum],0);
419
420 //delete inputFile;
421 // inputFile= NULL;
422 return tree;
423}
424
426
427 if ( m_fileNum >= int(m_fileNames.size())-1 ){
428 if(m_inputFiles[0]){
429 delete m_inputFiles[0];
430 m_inputFiles[0] = NULL;
431 }
432 return true;
433 }
434
435 (*m_branchesRead).Clear();
436 unsigned int treenr;
437 getTreeNr("Event",treenr);
438 if(m_inputFiles[treenr]){
439 delete m_inputFiles[treenr];
440 m_inputFiles[treenr] = NULL;
441 }
442 if(m_inputTrees[treenr]){
443 //delete m_inputTrees[treenr];
444 m_inputTrees[treenr] = NULL;
445 }
446 if(m_otherTrees[m_fileNum]) delete m_otherTrees[m_fileNum];
447 if(inputFiles[m_fileNum]) delete inputFiles[m_fileNum];
448
449 m_ENDFILE = true;
450 m_fileNum++;
451
452 m_entries=-1;
453 m_EOF = false;
454 return false;
455}
456
457/*
458bool RootInterface::checkEndOfTree(){
459 static int fileNum = 1;
460
461 if (m_EOF){
462 if (fileNmu >m_fileNames.size())
463 return true; //End of all Files;
464 else{
465 Tfile inputFile = TFile::Open(m_fileNames[fileNum].c_str(),"READ");
466 TTree *tree;
467 if(!inputFile->IsOpen()){
468 log<<MSG::ERROR<<"ROOT File" <<inputFile->GetName()<<"Coult not be opened for reading."<<endreq;
469 delete inputFile;
470 tree = Null; // The Root can not be opened
471 }
472 log<<MSG::INFO<<"RootIntrFace:;Opened File for input:"<<m_fileNames[fileNum].c_str()<<endreq;
473 tree =(TTree*)m_fileNames[fileNum]->Get("Dst");//the same tree name;
474 if (m_entries<=0){
475 m_entries=(Int_t)tree->GetEntries();
476 }
477 return false;
478 }
479 }
480 return false;
481}
482*/
483
484StatusCode RootInterface::setBranchAddress(const std::string treename, const std::string branchname, void *addr,int &branchnr)
485{
486 log << MSG::DEBUG <<"RootInterface::setbranch address, treename: "<<treename <<", branch "<<branchname <<endreq;
487
488 branchnr=-1;
489
490 TTree * tree;
491
492 if(m_fileNum != 0){
493
494 tree = getOtherTree(treename);
495 }else{
496 tree = getTree(treename);
497 }
498
499 //TTree * tree = getTree("Dst");
500 if (!tree) {
501 log << MSG::ERROR << "Could not find tree " << treename <<endreq;
502 log << MSG::ERROR << "terminate the process handly"<<endreq;
503 exit(1);
504 return StatusCode::FAILURE;
505 }
506 tree->SetMakeClass(1); // necessary for separate branch reading (segv otherwise)!
507
508 // log << MSG::INFO <<"ok!!!!!!!!!!!!!11"<<endreq;
509 TBranch *b = tree->GetBranch(branchname.c_str());
510 if (!b) {
511 //tree->Print();
512 log << MSG::DEBUG << "Could not find branch xx" << branchname <<"xx"<<endreq;
513 return StatusCode::FAILURE;
514 }
515 // log << MSG::INFO <<"ok!!!!!!!!!!!!!22"<<endreq;
516 // log << MSG::INFO <<"ok!!!!!!!!!!!!!22"<<(*addr)<<endreq;
517 b->SetAddress(addr);
518 // log << MSG::INFO <<"ok!!!!!!!!!!!!!33"<<endreq;
519 branchnr=m_branchesRead->GetEntries();
520 // log << MSG::INFO <<"ok!!!!!!!!!!!!!44"<<endreq;
521 TClonesArray &a = *m_branchesRead;
522 m_branchesRead->Expand(branchnr+1);
523 a[branchnr]=b;
524 return StatusCode::SUCCESS;
525}
526
527StatusCode RootInterface::getBranchEntry(int nr, int entry, void * addr, int& nb)
528{
529 log << MSG::DEBUG <<"RootInterface::getBranchEntry: "<<", branch nr "<<nr <<", entry "<<entry <<endreq;
530
531 if (nr <0) return StatusCode::FAILURE;
532 TBranch *branch=(TBranch *)m_branchesRead->At(nr);
533 if (!branch) {
534 log << MSG::ERROR << "Could not find branch " << nr <<endreq;
535 return StatusCode::FAILURE;
536 }
537
538 branch->SetAddress(addr);
539 nb=branch->GetEntry(entry);
540
541 if (nb<=0){
542 m_EOF=true;
543
544 }
545 return StatusCode::SUCCESS;
546}
547
548
549StatusCode RootInterface::getBranchEntry(int nr, int entry, int& nb)
550{
551 log << MSG::DEBUG <<"RootInterface::getBranchEntry: "<<", branch nr "<<nr <<", entry "<<entry <<endreq;
552
553 if (nr <0) return StatusCode::FAILURE;
554 TBranch *branch=(TBranch *)m_branchesRead->At(nr);
555 if (!branch) {
556 log << MSG::ERROR << "Could not find branch " << nr <<endreq;
557 return StatusCode::FAILURE;
558 }
559 nb=branch->GetEntry(entry);
560
561 if (nb<=0){
562 m_EOF=true;
563
564 }
565
566 return StatusCode::SUCCESS;
567}
568
569StatusCode RootInterface::getTreeNr(const std::string treename,unsigned int& treenr,bool doAdd) {
570 // look whether this tree has already got a number
571 // if not, add it
572 std::vector<std::string>::iterator where=std::find(m_treenames.begin(),m_treenames.end(),treename);
573 if (where == m_treenames.end()) {
574 if (doAdd) {
575 treenr=m_treenames.size();
576 m_treenames.push_back(treename);
577 m_inputFilenames.push_back("");
578 m_inputFiles.push_back(NULL);
579 m_inputTrees.push_back(NULL);
580 m_outputFilenames.push_back("");
581 m_outputFiles.push_back(NULL);
582 m_outputTrees.push_back(NULL);
583 m_splitModes.push_back(0);
584 m_bufSizes.push_back(0);
585 m_compressionLevels.push_back(0);
586 return StatusCode::SUCCESS;
587 }else {
588 log << MSG::ERROR << "Invalid tree name: " <<treename<< endreq;
589 return StatusCode::FAILURE;
590 }
591 }
592 treenr=where-m_treenames.begin();
593 return StatusCode::SUCCESS;
594}
595
597 // loop over all trees and fill them
598 StatusCode sc=StatusCode::FAILURE;
599 int nb;
600 std::vector<TTree *>::const_iterator trees;
601 for (trees=m_outputTrees.begin();trees<m_outputTrees.end();trees++) {
602 if ((*trees)==NULL) continue;
603 int treenr=(*trees)->GetUniqueID();
604 if(m_outputFiles[treenr]->IsZombie()||(!m_outputFiles[treenr]->IsOpen())){
605 std::cout<<"RootInterface ERROR::The ROOT File:"<<m_outputFilenames[treenr].c_str()<<"status is false"<<std::endl;
606 exit(1);
607 }
608 nb=(*trees)->Fill();
609 m_outputFiles[treenr] = (*trees)->GetCurrentFile();
610 log << MSG::DEBUG << "filled tree "<<(* trees)->GetName() <<" with "<<nb<<" bytes"<< endreq;
611 if(nb==-1){
612 log << MSG::FATAL << "Error in filling tree "<<(* trees)->GetName() <<" with "<<nb<<" bytes"<< endreq;
613 exit(1);
614 }
615 sc=StatusCode::SUCCESS;
616 }
617 return sc;
618}
619
620
621
622
623
624StatusCode RootInterface::f_addOutput(const std::string& treename,
625 const std::string& file,
626 int splitx, int bufsize, int compression){
627 log << MSG::INFO << "addOutput to single event" << endreq;
628 StatusCode status = StatusCode::FAILURE;
629 unsigned int treenr;
630
631 status = f_getTreeNr(treename, treenr, true);
632 m_single_compressionLevels[treenr] = compression;
633 m_single_outputFileNames[treenr] = file;
634 m_single_outputFiles[treenr] = NULL;
635 m_single_outputTrees[treenr] = NULL;
636 m_single_splitModes[treenr] = splitx;
637 m_single_bufSizes[treenr] = bufsize;
638
639 std::cout << "finish f_addOutput to single event" << std::endl;
640 return status;
641}
642
643StatusCode RootInterface::f_createTree(unsigned int treenr,
644 const std::string treename){
645 log << MSG::INFO << "f_createTree()" << endreq;
646
647 TDirectory *saveDir = gDirectory;
648
649 m_single_outputFiles[treenr] =
650 new TFile(m_single_outputFileNames[treenr].c_str(), "RECREATE");
651 if ( !m_single_outputFiles[treenr]->IsOpen()){
652 log << MSG::ERROR << "ROOT share file: "
653 << m_single_outputFileNames[treenr]
654 << " could not be opened for writing"
655 << endreq;
656 return StatusCode::FAILURE;
657 }
658 log << MSG::INFO << "f_createTree()::open share file for writing: "
659 << m_single_outputFileNames[treenr] << endreq;
660
661 m_single_outputFiles[treenr]->cd();
662 m_single_outputFiles[treenr]->SetCompressionLevel(m_single_compressionLevels[treenr]);
663
664 std::string title = treename + " for share";
665 m_single_outputTrees[treenr] = new TTree(treename.c_str(), title.c_str());
666 saveDir->cd();
667
668 return StatusCode::SUCCESS;
669}
670
671StatusCode RootInterface::f_createBranch(const std::string& treename,
672 const std::string& branchname,
673 const char *classname,
674 void *addr, int & branchnr){
675 log << MSG::INFO << "f_craeteBranch() create branch, tree name:"
676 << treename << ", branch name:" << branchname << endreq;
677
678 TBranch *branch;
679 unsigned int treenr;
680 StatusCode status = f_getTreeNr(treename, treenr);
681 if ( !status.isSuccess()) return status;
682
683 if ( !m_single_outputTrees[treenr])
684 status = this->f_createTree(treenr, treename);
685 if ( !status.isSuccess()) return status;
686
687 TTree* tree = m_single_outputTrees[treenr];
688 tree->SetUniqueID(treenr);
689
690 branch = tree->Branch(branchname.c_str(),
691 classname,
692 addr,
693 m_single_bufSizes[treenr],
694 m_single_splitModes[treenr]);
695
696}
697
698StatusCode RootInterface::f_getTreeNr(const std::string treename,
699 unsigned int& treenr,bool doAdd){
700
701 std::vector<std::string>::iterator where =
702 std::find(m_single_treenames.begin(), m_single_treenames.end(), treename);
703
704 if ( where == m_single_treenames.end()){
705 if ( doAdd){
706 treenr = m_single_treenames.size();
707 m_single_treenames.push_back(treename);
708
709 m_single_outputFileNames.push_back("");
710 m_single_outputFiles.push_back(NULL);
711 m_single_outputTrees.push_back(NULL);
712 m_single_splitModes.push_back(0);
713 m_single_bufSizes.push_back(0);
714 m_single_compressionLevels.push_back(0);
715
716 return StatusCode::SUCCESS;
717 }
718 else {
719 log << MSG::ERROR << "Invalid share tree name: "
720 << treename << endreq;
721 return StatusCode::FAILURE;
722 }
723 }
724 treenr = where - m_single_treenames.begin();
725 return StatusCode::SUCCESS;
726}
727
729 StatusCode status = StatusCode::FAILURE;
730 int byte;
731
732 std::vector<TTree *>::const_iterator tree;
733 for ( tree = m_single_outputTrees.begin(); tree < m_single_outputTrees.end(); tree++){
734 if ( (*tree) == NULL) continue;
735 byte = (*tree)->Fill();
736 (*tree)->Print();
737 log << MSG::INFO << "f_fillTrees() filled tree " << (*tree)->GetName()
738 << " with " << byte << " bytes!" << endreq;
739 status = StatusCode::SUCCESS;
740 }
741
742 return status;
743}
744
746 log << MSG::INFO << "f_finalize() in RootInterface" << endreq;
747
748 std::vector<TTree *>::const_iterator tree;
749 for ( tree = m_single_outputTrees.begin(); tree < m_single_outputTrees.end(); tree++){
750 if ( *tree){
751 unsigned int treenr = (*tree)->GetUniqueID();
752 log << MSG::INFO << "tree id: " << treenr << endreq;
753 if ( m_single_outputFiles[treenr] ){
754 if ( !m_single_outputFiles[treenr]->IsOpen()){
755 log << MSG::ERROR << "f_finalize could not open share file for writing"
756 << endreq;
757 return StatusCode::FAILURE;
758 }
759 else {
760 log << MSG::INFO << "Closing file:" << treenr
761 << ", tree:" << (*tree)->GetName() << endreq;
762
763 TDirectory *saveDir = gDirectory;
764 m_single_outputFiles[treenr]->cd();
765 log <<MSG::INFO << "WREITE TO FILE BYTES: "
766 << m_single_outputFiles[treenr]->Write()
767 << endreq;
768 m_single_outputFiles[treenr]->Close();
769 saveDir->cd();
770 }
771 }
772 }
773 }
774 return StatusCode::SUCCESS;
775}
string getDecayOptions()
Definition: DataInfoSvc.h:25
void setTotEvtNo(std::vector< int > i)
Definition: DataInfoSvc.h:29
std::vector< int > getTotEvtNo()
Definition: DataInfoSvc.h:26
virtual std::vector< int > getTotEvtNo()
RootInterface(MsgStream log)
virtual StatusCode f_createBranch(const std::string &treename, const std::string &branchname, const char *classname, void *addr, int &branchnr)
virtual StatusCode f_getTreeNr(const std::string treename, unsigned int &treenr, bool doAdd=false)
virtual StatusCode f_finalize()
virtual std::string getJobOptions()
virtual StatusCode getBranchEntry(int nr, int entry, int &nb)
get entry from this branch
virtual StatusCode f_fillTrees()
virtual StatusCode finalize()
virtual bool checkEndOfTree()
check if all the files is over 2005-11-28
virtual StatusCode createBranch(const std::string &tree, const std::string &branch, const char *classname, void *addr, int &branchnr)
create a branch in this tree
virtual void printJobInfo(TFile *file, int level)
virtual StatusCode setBranchAddress(const std::string treename, const std::string branchname, void *addr, int &nb)
set branch address
virtual StatusCode addInput(const std::string &treename, const std::string &file)
add input tree to the list
virtual StatusCode f_addOutput(const std::string &treename, const std::string &file, int splitx=1, int bufsize=64000, int compression=1)
virtual ~RootInterface()
static RootInterface * Instance(MsgStream log)
singleton behaviour
virtual StatusCode addOutput(const std::string &treename, const std::string &file, int splitx, int bufsize, int compression)
add output tree to the list
virtual StatusCode fillTrees()
fill in all trees
virtual StatusCode f_createTree(unsigned int treenr, const std::string treename)
virtual std::string getDecayOptions()
void setBossVer(string ver)
Definition: TJobInfo.h:21
void setTotEvtNo(std::vector< int > i)
Definition: TJobInfo.h:25
string getDecayOptions() const
Definition: TJobInfo.h:18
void setDecayOptions(string opt)
Definition: TJobInfo.h:24
vector< string > getJobOptions() const
Definition: TJobInfo.h:17
std::vector< int > getTotEvtNo() const
Definition: TJobInfo.h:19
string getBossVer() const
Definition: TJobInfo.h:16
void setJobOptions(vector< string > opt)
Definition: TJobInfo.h:23
char * c_str(Index i)
Definition: EvtCyclic3.cc:252