The Gaudi Framework  master (ff829712)
Loading...
Searching...
No Matches
Gaudi::RootDatabaseMerger Class Reference
Collaboration diagram for Gaudi::RootDatabaseMerger:

Classes

union  uuid_data
 

Public Member Functions

 RootDatabaseMerger ()
 Standard constructor.
 
virtual ~RootDatabaseMerger ()
 Default destructor.
 
bool exists (const std::string &fid) const
 Check if a database exists.
 
MergeStatus create (const std::string &fid)
 Create new output file.
 
MergeStatus attach (const std::string &fid)
 Attach to existing output file for further merging.
 
MergeStatus close ()
 Close output file.
 
MergeStatus merge (const std::string &fid)
 Merge new input to existing output.
 
void dumpSections ()
 Dump collected database sections.
 
MergeStatus saveSections ()
 Save new sections to the output file.
 
MergeStatus createFID ()
 Create and add new FID to the newly merged file.
 
MergeStatus copyAllTrees (TFile *source)
 Copy all data trees from the input file to the output file.
 
MergeStatus copyTree (TFile *source, const std::string &name)
 Copy one single tree from the input file to the output file.
 
MergeStatus copyBranch (TTree *src_tree, TTree *out_tree, const std::string &name)
 Copy single reference branch.
 
MergeStatus copyRefs (TFile *source, const std::string &name)
 Copy one single tree from the input file to the output file.
 
MergeStatus addSections (TTree *in, TTree *out)
 Add section information for the next merge step.
 

Private Types

typedef std::vector< ContainerSectionContainerSections
 
typedef std::map< std::string, ContainerSectionsDatabaseSections
 

Private Attributes

DatabaseSections m_sections
 
std::unique_ptr< TFile > m_output
 
bool m_treeSections = false
 

Detailed Description

Author
M.Frank
Version
1.0

Definition at line 50 of file merge.C.

Member Typedef Documentation

◆ ContainerSections

Definition at line 52 of file merge.C.

◆ DatabaseSections

typedef std::map<std::string, ContainerSections> Gaudi::RootDatabaseMerger::DatabaseSections
private

Definition at line 53 of file merge.C.

Constructor & Destructor Documentation

◆ RootDatabaseMerger()

RootDatabaseMerger::RootDatabaseMerger ( )
default

Standard constructor.

◆ ~RootDatabaseMerger()

RootDatabaseMerger::~RootDatabaseMerger ( )
virtual

Default destructor.

Definition at line 111 of file merge.C.

111{ close(); }
MergeStatus close()
Close output file.
Definition merge.C:194

Member Function Documentation

◆ addSections()

MergeStatus RootDatabaseMerger::addSections ( TTree * in,
TTree * out )

Add section information for the next merge step.

Definition at line 289 of file merge.C.

289 {
290 if ( m_treeSections ) {
291 ContainerSection s;
292 s.start = (int)( out ? out->GetEntries() : 0 );
293 s.length = (int)in->GetEntries();
294 m_sections[in->GetName()].push_back( s );
295 return MERGE_SUCCESS;
296 }
297
298 TObjArray* a_in = in->GetListOfBranches();
299 for ( int i = 0, n = a_in->GetLast(); i < n; ++i ) {
300 TBranch* b_in = (TBranch*)a_in->At( i );
301 TBranch* b_out = out ? out->GetBranch( b_in->GetName() ) : 0;
302 if ( !out || b_out ) {
303 ContainerSection s;
304 s.start = (int)( b_out ? b_out->GetEntries() : 0 );
305 s.length = (int)b_in->GetEntries();
306 m_sections[b_in->GetName()].push_back( s );
307 continue;
308 }
309 ::printf( "+++ Cannot merge incompatible branches:%s.\n", b_in->GetName() );
310 return MERGE_ERROR;
311 }
312 return MERGE_SUCCESS;
313}
DatabaseSections m_sections
Definition merge.C:55

◆ attach()

MergeStatus RootDatabaseMerger::attach ( const std::string & fid)

Attach to existing output file for further merging.

Definition at line 121 of file merge.C.

121 {
122 const char* fid = file_id.c_str();
123 if ( m_output ) {
124 ::printf( "+++ Another output file %s is already open. Request denied.\n", m_output->GetName() );
125 return MERGE_ERROR;
126 } else if ( !exists( file_id ) ) {
127 ::printf( "+++ Cannot attach output file %s --- file does not exist.\n", fid );
128 return MERGE_ERROR;
129 }
130 m_output.reset( TFile::Open( fid, "UPDATE" ) );
131 if ( m_output && !m_output->IsZombie() ) {
132 // if ( s_dbg ) ::printf("+++ Update output file %s.\n",fid);
133 // if ( s_dbg ) ::printf("+++ Update output file.\n");
134 return MERGE_SUCCESS;
135 }
136 ::printf( "+++ Failed to open new output file %s.\n", fid );
137 return MERGE_ERROR;
138}
std::unique_ptr< TFile > m_output
Definition merge.C:56
bool exists(const std::string &fid) const
Check if a database exists.
Definition merge.C:114

◆ close()

MergeStatus RootDatabaseMerger::close ( )

Close output file.

Definition at line 194 of file merge.C.

194 {
195 if ( m_output ) {
196 // if ( s_dbg ) ::printf("+++ Closing merge file: %s.\n",m_output->GetName());
197 if ( s_dbg ) ::printf( "+++ Closing merge file.\n" );
198 m_output->Write();
199 m_output->Purge();
200 if ( s_dbg ) m_output->ls();
201 m_output->Close();
202 m_output.reset();
203 }
204 return MERGE_SUCCESS;
205}

◆ copyAllTrees()

MergeStatus RootDatabaseMerger::copyAllTrees ( TFile * source)

Copy all data trees from the input file to the output file.

Definition at line 333 of file merge.C.

333 {
334 TIter nextkey( source->GetListOfKeys() );
335 // m_treeSections = true;
336 for ( TKey* key = (TKey*)nextkey(); key; key = (TKey*)nextkey() ) {
337 const char* classname = key->GetClassName();
338 TClass* cl = gROOT->GetClass( classname );
339 if ( !cl ) continue;
340 if ( cl->InheritsFrom( "TTree" ) ) {
341 string name = key->GetName();
342 if ( name == "Refs" ) continue;
343 m_treeSections = 0 == ::strncmp( key->GetName(), "<local>_", 7 );
344 printf( "+++ Copy Tree:%s %d\n", name.c_str(), int( m_treeSections ) );
345 if ( copyTree( source, name ) != MERGE_SUCCESS ) {
346 m_treeSections = false;
347 return MERGE_ERROR;
348 }
349 }
350 }
351 m_treeSections = false;
352 return MERGE_SUCCESS;
353}
MergeStatus copyTree(TFile *source, const std::string &name)
Copy one single tree from the input file to the output file.
Definition merge.C:356

◆ copyBranch()

MergeStatus RootDatabaseMerger::copyBranch ( TTree * src_tree,
TTree * out_tree,
const std::string & name )

Copy single reference branch.

Definition at line 316 of file merge.C.

316 {
317 char text[4096];
318 TBranch* s = src_tree->GetBranch( name.c_str() );
319 TBranch* o = out_tree->GetBranch( name.c_str() );
320 if ( s && o ) {
321 s->SetAddress( text );
322 o->SetAddress( text );
323 for ( Long64_t i = 0, n = s->GetEntries(); i < n; ++i ) {
324 s->GetEntry( i );
325 o->Fill();
326 }
327 return MERGE_SUCCESS;
328 }
329 return MERGE_ERROR;
330}

◆ copyRefs()

MergeStatus RootDatabaseMerger::copyRefs ( TFile * source,
const std::string & name )

Copy one single tree from the input file to the output file.

Copy refs of one single tree from the input file to the output file.

Definition at line 400 of file merge.C.

400 {
401 TTree* src_tree = (TTree*)source->Get( name.c_str() );
402 if ( src_tree ) {
403 TTree* out_tree = (TTree*)m_output->Get( name.c_str() );
404 if ( out_tree ) {
405 addSections( src_tree, out_tree );
406 copyBranch( src_tree, out_tree, "Links" );
407 copyBranch( src_tree, out_tree, "Params" );
408 copyBranch( src_tree, out_tree, "Containers" );
409 copyBranch( src_tree, out_tree, "Databases" );
410 out_tree->Write();
411 return MERGE_SUCCESS;
412 }
413 }
414 return MERGE_ERROR;
415}
MergeStatus addSections(TTree *in, TTree *out)
Add section information for the next merge step.
Definition merge.C:289
MergeStatus copyBranch(TTree *src_tree, TTree *out_tree, const std::string &name)
Copy single reference branch.
Definition merge.C:316

◆ copyTree()

MergeStatus RootDatabaseMerger::copyTree ( TFile * source,
const std::string & name )

Copy one single tree from the input file to the output file.

Definition at line 356 of file merge.C.

356 {
357 TTree* src_tree = (TTree*)source->Get( name.c_str() );
358 if ( src_tree ) {
359 Long64_t src_entries = src_tree->GetEntries();
360 TTree* out_tree = (TTree*)m_output->Get( name.c_str() );
361 m_output->cd();
362 if ( 0 == src_tree->GetEntries() ) { src_tree->SetEntries( 1 ); }
363 addSections( src_tree, out_tree );
364 if ( !out_tree ) {
365 out_tree = src_tree->CloneTree( -1, "fast" );
366 if ( s_dbg ) ::printf( "+++ Created new Tree %s.\n", out_tree->GetName() );
367 out_tree->Write();
368 delete out_tree;
369 return MERGE_SUCCESS;
370 }
371 m_output->GetObject( name.c_str(), out_tree );
372 TTreeCloner cloner( src_tree, out_tree, "fast" );
373 if ( cloner.IsValid() ) {
374 Long64_t out_entries = out_tree->GetEntries();
375 out_tree->SetEntries( out_entries + src_entries );
376 Bool_t res = cloner.Exec();
377 if ( s_dbg ) ::printf( "+++ Merged tree: %s res=%d\n", out_tree->GetName(), res );
378 out_tree->Write();
379 delete out_tree;
380 return MERGE_SUCCESS;
381 } else {
382 // Fast cloning is not possible for this input TTree.
383 // ... see TTree::CloneTree for example of recovery code ...
384 ::printf( "+++ Got a tree where fast cloning is not possible -- operation failed.\n" );
385 return MERGE_ERROR;
386 }
387#if 0
388 Long64_t nb = out_tree->CopyEntries(src_tree,-1,"fast");
389 Long64_t out_entries = out_tree->GetEntries();
390 out_tree->SetEntries(out_entries+src_entries);
391 out_tree->Write();
392 if ( s_dbg ) ::printf("+++ Merged tree: %s res=%lld\n",out_tree->GetName(),nb);
393 return MERGE_SUCCESS;
394#endif
395 }
396 return MERGE_ERROR;
397}

◆ create()

MergeStatus RootDatabaseMerger::create ( const std::string & fid)

Create new output file.

Definition at line 141 of file merge.C.

141 {
142 if ( m_output ) {
143 ::printf( "+++ Another output file %s is already open. Request denied.\n", m_output->GetName() );
144 return MERGE_ERROR;
145 } else if ( exists( fid ) ) {
146 ::printf( "+++ Cannot create output file %s --- file already exists.\n", fid.c_str() );
147 return MERGE_ERROR;
148 }
149 m_output.reset( TFile::Open( fid.c_str(), "RECREATE" ) );
150 if ( m_output && !m_output->IsZombie() ) {
151 TTree* t1 = new TTree( "Sections", "Root Section data" );
152 TTree* t2 = new TTree( "Refs", "Root Section data" );
153 if ( t1 && t2 ) {
154 t1->Branch( "Sections", 0, "Sections/C" );
155 t2->Branch( "Links", 0, "Links/C" );
156 t2->Branch( "Params", 0, "Params/C" );
157 t2->Branch( "Databases", 0, "Databases/C" );
158 t2->Branch( "Containers", 0, "Containers/C" );
159 if ( s_dbg ) ::printf( "+++ Opened new output file %s.\n", fid.c_str() );
160 return MERGE_SUCCESS;
161 }
162 }
163 ::printf( "+++ Failed to open new output file %s.\n", fid.c_str() );
164 return MERGE_ERROR;
165}
TemplatedAlg< int, std::vector< std::string > > t1
TemplatedAlg< double, bool > t2

◆ createFID()

MergeStatus RootDatabaseMerger::createFID ( )

Create and add new FID to the newly merged file.

Close output file.

Definition at line 168 of file merge.C.

168 {
169 static const char* fmt = "FID=%08lX-%04hX-%04hX-%02hhX%02hhX-%02hhX%02hhX%02hhX%02hhX%02hhX%02hhX";
170 if ( m_output ) {
171 TTree* t = (TTree*)m_output->Get( "Refs" );
172 if ( t ) {
173 uuid_data d;
174 char text[256];
175 TUUID uuid;
176 TBranch* b = t->GetBranch( "Params" );
177 if ( b ) {
178 uuid.GetUUID( d.buf );
179 sprintf( text, fmt, d.ibuf[0], d.sbuf[2], d.sbuf[3], d.buf[8], d.buf[9], d.buf[10], d.buf[11], d.buf[12],
180 d.buf[13], d.buf[14], d.buf[15] );
181 b->SetAddress( text );
182 b->Fill();
183 t->Write();
184 if ( s_dbg ) ::printf( "+++ Added new GUID %s to merge file.\n", text );
185 return MERGE_SUCCESS;
186 }
187 }
188 }
189 ::printf( "+++ Failed to add new GUID to merge file.\n" );
190 return MERGE_ERROR;
191}

◆ dumpSections()

void RootDatabaseMerger::dumpSections ( )

Dump collected database sections.

Definition at line 251 of file merge.C.

251 {
252 for ( const auto& i : m_sections ) {
253 int cnt = 0;
254 string prefix = i.first;
255 for ( const auto& j : i.second ) {
256 char text[1024];
257 ::sprintf( text, "['%s'][%d]", prefix.c_str(), cnt++ );
258 if ( s_dbg ) {
259 ::printf( "+++ section %-55s Start:%8d ... %8d [%d entries]\n", text, j.start, j.start + j.length, j.length );
260 }
261 }
262 }
263}
str prefix
Definition gaudirun.py:361

◆ exists()

bool RootDatabaseMerger::exists ( const std::string & fid) const

Check if a database exists.

Definition at line 114 of file merge.C.

114 {
115 Bool_t result = gSystem->AccessPathName( fid.c_str(), kFileExists );
116 // if ( s_dbg ) ::printf("File %s %s!\n",fid.c_str(),result == kFALSE ? "EXISTS" : "DOES NOT EXIST");
117 return result == kFALSE;
118}

◆ merge()

MergeStatus RootDatabaseMerger::merge ( const std::string & fid)

Merge new input to existing output.

Definition at line 266 of file merge.C.

266 {
267 if ( m_output ) {
268 std::unique_ptr<TFile> source{ TFile::Open( fid.c_str() ) };
269 if ( source && !source->IsZombie() ) {
270 size_t idx = fid.rfind( '/' );
271 ::printf( "+++ Start merging input file:%s\n",
272 idx != string::npos ? fid.substr( idx + 1 ).c_str() : fid.c_str() );
273 if ( copyAllTrees( source.get() ) == MERGE_SUCCESS ) {
274 if ( copyRefs( source.get(), "Refs" ) == MERGE_SUCCESS ) {
275 source->Close();
276 return MERGE_SUCCESS;
277 }
278 }
279 }
280 ::printf( "+++ Cannot open input file:%s\n", fid.c_str() );
281 m_output->cd();
282 return MERGE_ERROR;
283 }
284 ::printf( "+++ No valid output file present. Merge request refused for fid:%s.\n", fid.c_str() );
285 return MERGE_ERROR;
286}
MergeStatus copyRefs(TFile *source, const std::string &name)
Copy one single tree from the input file to the output file.
Definition merge.C:400
MergeStatus copyAllTrees(TFile *source)
Copy all data trees from the input file to the output file.
Definition merge.C:333

◆ saveSections()

MergeStatus RootDatabaseMerger::saveSections ( )

Save new sections to the output file.

Definition at line 208 of file merge.C.

208 {
209 if ( m_output ) {
210 int nb, total = 0, nbytes = 0;
211 char text[1024];
212 TTree* t = (TTree*)m_output->Get( "Sections" );
213 if ( t ) {
214 TBranch* b = t->GetBranch( "Sections" );
215 if ( b ) {
216 b->SetAddress( text );
217 for ( const auto& i : m_sections ) {
218 string cont = i.first;
219 for ( const auto& j : i.second ) {
220 ::sprintf( text, "[CNT=%s][START=%d][LEN=%d]", cont.c_str(), j.start, j.length );
221 nb = t->Fill();
222 ++total;
223 if ( nb > 0 )
224 nbytes += nb;
225 else
226 ::printf( "+++ Failed to update Sections tree with new entries. [WRITE_ERROR]\n" );
227 }
228 }
229 ::sprintf( text, "[END-OF-SECTION]" );
230 nb = t->Fill();
231 ++total;
232 if ( nb > 0 )
233 nbytes += nb;
234 else
235 ::printf( "+++ Failed to update Sections branch with new entries. [WRITE_ERROR]\n" );
236 t->Write();
237 if ( s_dbg ) ::printf( "+++ Added %d Sections entries with %d bytes in total.\n", total, nbytes );
238 return MERGE_SUCCESS;
239 }
240 ::printf( "+++ Failed to update Sections tree with new entries. [NO_OUTPUT_BRANCH]\n" );
241 return MERGE_ERROR;
242 }
243 ::printf( "+++ Failed to update Sections tree with new entries. [NO_OUTPUT_TREE]\n" );
244 return MERGE_ERROR;
245 }
246 ::printf( "+++ Failed to update Sections tree with new entries. [NO_OUTPUT_FILE]\n" );
247 return MERGE_ERROR;
248}

Member Data Documentation

◆ m_output

std::unique_ptr<TFile> Gaudi::RootDatabaseMerger::m_output
private

Definition at line 56 of file merge.C.

◆ m_sections

DatabaseSections Gaudi::RootDatabaseMerger::m_sections
private

Definition at line 55 of file merge.C.

◆ m_treeSections

bool Gaudi::RootDatabaseMerger::m_treeSections = false
private

Definition at line 57 of file merge.C.


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