]> Creatis software - clitk.git/blob - tests_dav/GateMergeManager.cc
switched to gjm merger for pet
[clitk.git] / tests_dav / GateMergeManager.cc
1 /*----------------------
2    GATE version name: gate_v...
3
4    Copyright (C): OpenGATE Collaboration
5
6 This software is distributed under the terms
7 of the GNU Lesser General  Public Licence (LGPL)
8 See GATE/LICENSE.txt for further details
9 ----------------------*/
10
11
12 #include <TROOT.h>
13 #include <TFile.h>
14 #include <TChain.h>
15 #include <TTree.h>
16 #include <TBranch.h>
17 #include <TIterator.h>
18 #include <TObject.h>
19 #include <TKey.h>
20 #include <TH1.h>
21 #include <TH2.h>
22
23 #include <iostream>
24 #include <fstream>
25 #include <string>
26 #include <glob.h>
27 #include <ctype.h>
28 #include <sstream>
29 #include <vector>
30 #include <cstdlib>
31 #include <cmath>
32
33 #include "GateMergeManager.hh"
34
35 using namespace std;
36
37 void GateMergeManager::StartMerging(string splitfileName){
38
39   // get the files to merge
40   ReadSplitFile(splitfileName);
41   //do the merging
42   if (m_fastMerge==true) FastMergeRoot();
43   else MergeRoot();
44
45   //if we are here the merging has been successful
46   //we mark the directory as ready for cleanup
47   //string dir   = splitfileName.substr(0,splitfileName.rfind("/"));
48   //string ready = "touch "+dir+"/ready_for_delete";
49   //const int res=system(ready.c_str());
50   //if(res) cout<<"Strange?? Can't mark "<<dir<<" as done!"<<endl;
51 };
52
53 // to process the splitfile
54 void GateMergeManager::ReadSplitFile(string splitfileName){
55
56   ifstream splitfile;
57   splitfile.open(splitfileName.c_str());
58   if(!splitfile){
59      cout<<"Can't open split file: "<<splitfileName<<endl;
60      exit(0);
61   };
62
63   //avoid multiple calls
64   if(m_vRootFileNames.size()>0) return;
65
66   stringstream ssnfiles;
67   char cline[512];
68   while(splitfile){               //we allow for rubbish before "Number of files:"
69      splitfile.getline(cline,512);
70      if(!strncmp(cline,"Number of files:",16)){
71          ssnfiles.str(cline+16);
72          ssnfiles>>m_Nfiles;      // the number of files to merge
73          break;
74      }
75   }
76
77   // now we look for the file names
78   int iRoot=0;
79   while(splitfile){
80      splitfile.getline(cline,512);
81      // check for root
82      // input files
83      if(!strncmp(cline,"Root filename:",14)){
84         m_vRootFileNames.push_back(strtok(cline+15," "));
85         m_vRootFileNames[iRoot]+=".root";
86         if(m_verboseLevel>2) cout<<"Root input file name: "<<m_vRootFileNames[iRoot]<<endl;
87         iRoot++;
88      }
89      // output file
90      else if(!strncmp(cline,"Original Root filename:",23)){
91         m_RootTargetName=strtok(cline+23," ");
92         m_RootTargetName+=".root";
93         if(m_outDir!=""){// output directory option used
94            //replace path with outDir
95            size_t pos=m_RootTargetName.rfind('/',m_RootTargetName.length());
96            if(pos==string::npos) pos=-1;
97            m_RootTargetName=m_outDir+m_RootTargetName.substr(pos+1);
98         }
99         if(m_verboseLevel>2) cout<<"Root output file name: "<<m_RootTargetName<<endl;
100      }
101   }
102
103   // check if number of root files correct
104   if(iRoot && iRoot!=m_Nfiles) {
105        cout<<"Inconsistent number of root file entries in split file!"<<endl;
106        exit(0);
107   }
108 }
109
110 /************************************************************************************/
111 void GateMergeManager::FastMergeRoot()
112 {
113    //check target name
114    if(m_RootTargetName=="") {
115       cout<<"No root output file name given in split file"<<endl;
116       exit(0);
117    }
118
119    // number of files to merge
120    int nfiles=m_vRootFileNames.size();
121    if(nfiles<1){
122       cout<<"No files to merge!"<<endl;
123       exit(0);
124    }
125    //we try to recover the last_event_ID in all root files
126    for(int i=0;i<nfiles;i++) {
127       filearr.push_back(TFile::Open(m_vRootFileNames[i].c_str(),"OLD"));
128       if(filearr[i]==NULL){
129          cout<<"Not a readable file "<<m_vRootFileNames[i]<<" - exit!"<<endl;   
130          exit(0);
131       }
132    }  
133    m_lastEvents.resize(nfiles+1);
134    m_lastEvents[0]=0;
135    bool flgLstEvntID(false);
136
137    // latest_event_ID histogram
138    flgLstEvntID=true;
139    TH1* h=NULL;
140    for(int i=0;i<nfiles;i++){
141      h = (TH1*)filearr[i]->Get("latest_event_ID");
142      m_lastEvents[i+1]=int(h->GetMean());
143      }
144
145    if(!flgLstEvntID) {
146       cout<<"FastMerge: No latest_event_ID histogram - exit!"<<endl;
147       exit(0);
148    }
149
150  /*for (int i=0;i<nfiles;i++)
151   {
152    cout<<"Last ID is: "<<m_lastEvents[i+1]<<" File " <<i+1<<endl;
153   }*/
154
155    //find all trees
156    vector<string> treeNames;
157    TFile* node = TFile::Open(m_vRootFileNames[0].c_str(),"OLD");
158    TIter nextkey(node->GetListOfKeys());
159    TKey *key=0;
160    while ((key = (TKey*)nextkey()))
161    {
162      node->cd();
163      TObject* obj = (TObject*)key->ReadObj();
164      if(  obj->IsA()->InheritsFrom("TTree"))
165      {
166        //no doubles
167        bool singleName=true;
168        for(unsigned int i=0;i<treeNames.size();i++)
169        {
170          if(treeNames[i].compare(obj->GetName())==0)
171          {
172            singleName=false;
173            break;
174          }
175        }
176        if(singleName) treeNames.push_back(obj->GetName());
177      }
178    }
179    //take care of all trees
180    for(unsigned int i=0;i<treeNames.size();i++) 
181    {
182      if(!MergeTree(treeNames[i])) if(m_verboseLevel>1) cout<<"Problem with merging "<<treeNames[i]<<endl; 
183    }
184 }
185 /************************************************************************************/
186 void GateMergeManager::MergeRoot(){
187
188    //check target name
189    if(m_RootTargetName=="") {
190       cout<<"No root output file name given in split file"<<endl;
191       exit(0);
192    }
193
194    // number of files to merge
195    int nfiles=m_vRootFileNames.size();
196
197    // the root output file
198    if(m_forced) {
199       m_RootTarget = TFile::Open(m_RootTargetName.c_str(),"RECREATE");
200    } else {
201       m_RootTarget = TFile::Open(m_RootTargetName.c_str(),"NEW");
202       if(!m_RootTarget){
203           cout<<"The root ouput file already exists! Try -f to overwrite and erase the file."<<endl;
204           exit(0);
205       }
206    }
207
208    if(nfiles<1){
209       cout<<"No files to merge!"<<endl;
210       exit(0);
211    }
212
213    TFile* filearr[nfiles];
214    for(int i=1;i<nfiles;i++) {
215       filearr[i] = TFile::Open(m_vRootFileNames[i].c_str(),"OLD");
216       if(filearr[i]==NULL){
217          cout<<"Not a readable file "<<m_vRootFileNames[i]<<" - exit!"<<endl;   
218          exit(0);
219       }
220    }
221
222 // first we copy all histos (only top directory) 
223 // and look for the latestEventID
224 // and find out which  trees/ntuples  exist
225    m_lastEvents.resize(nfiles+1);
226        m_lastEvents[0]=0;
227
228    bool flgLstEvntID(false);
229
230    if(m_verboseLevel>0) {
231       cout<<"Combining: ";
232       for(int i=0;i<nfiles;i++) cout<<m_vRootFileNames[i]<<" ";
233       cout<<"-> "<<m_RootTarget->GetName()<<endl;
234    }
235
236    vector<string> treeNames;
237    TFile* node = TFile::Open(m_vRootFileNames[0].c_str(),"OLD");
238    if(node==NULL){
239        cout<<"Not a readable file "<<m_vRootFileNames[0]<<" - exit!"<<endl;   
240        exit(0);
241    }
242         TIter nextkey(node->GetListOfKeys());
243         TKey *key=0;
244         while ((key = (TKey*)nextkey())) {
245           node->cd();
246           TObject* obj = (TObject*)key->ReadObj();
247           if(  obj->IsA()->InheritsFrom("TH1")){
248              m_RootTarget->cd();
249              TH1 *h1 = (TH1 *)obj->Clone();
250              if(strcmp(h1->GetName(),"latest_event_ID")){
251                      //any other histogram
252                      for(int i=1;i<nfiles;i++){
253                           TH1 *h2 = (TH1*)filearr[i]->Get( h1->GetName() );
254                           h1->Add( h2 ); 
255                      }
256                      h1->Write();
257              } else {
258                           // latest_event_ID histogram
259                           flgLstEvntID=true;
260                           m_lastEvents[1]=int(h1->GetMean());
261                           TH1* h=NULL;
262                           for(int i=1;i<nfiles;i++){
263                                   h = (TH1*)filearr[i]->Get("latest_event_ID");
264                                   m_lastEvents[i+1]=int(h->GetMean());
265                           }
266                           if(h!=NULL) h->Write(); // we keep the highest latest_event_ID to allow for successive merging
267                           else h1->Write();
268               }
269          }
270          if(  obj->IsA()->InheritsFrom("TH2")){
271            m_RootTarget->cd();
272            TH2 *h1 = (TH2 *)obj->Clone();
273                         for(int i=1;i<nfiles;i++){
274                                 TH2 *h2 = (TH2*)filearr[i]->Get( h1->GetName() );
275                                 h1->Add( h2 ); 
276                         }
277                         h1->Write();
278                 }
279          if(  obj->IsA()->InheritsFrom("TTree")){
280            //no doubles
281            bool singleName=true;
282            for(unsigned int i=0;i<treeNames.size();i++) {
283               if(treeNames[i].compare(obj->GetName())==0){
284                  singleName=false;
285                  break;
286               }
287            }
288            if(singleName) treeNames.push_back(obj->GetName());
289          }
290    }
291    if(!flgLstEvntID) {
292       cout<<"Merge: No latest_event_ID histogram - exit!"<<endl;
293       exit(0);
294    }
295
296    //now we take care of the trees
297    for(unsigned int i=0;i<treeNames.size();i++) 
298      if(!MergeTree(treeNames[i])) if(m_verboseLevel>1) cout<<"Problem with merging "<<treeNames[i]<<endl; 
299
300    // everything is done
301 };
302
303 /*******************************************************************************************/
304 // cleanup after merging
305 void GateMergeManager::StartCleaning(string splitfileName,bool test){
306
307   // get the files to erase
308   ReadSplitFile(splitfileName);
309
310   string dir=splitfileName.substr(0,splitfileName.rfind("/"));
311
312   // test for the mark: ready_for_delete 
313   string touched=dir+"/ready_for_delete";
314   ifstream ready;
315   ready.open(touched.c_str());
316   if (!ready) {
317      cout<<"Cannot do the cleanup - directory "<<dir<<" not marked as ready!"<<endl;
318      if(!test) exit(0);
319   }
320
321   // print info
322   if(test)                    cout<<"I would like to erase:"<<endl;
323   if(!test&&m_verboseLevel>1) cout<<"Going to erase the following files:"<<endl;
324   if(m_verboseLevel>1||test){
325      cout<<" --> "<<dir<<endl;
326      // root
327      for(unsigned int i=0;i<m_vRootFileNames.size();i++){
328         cout<<" --> "<<m_vRootFileNames[i]<<endl;
329      }
330   }
331
332   // erase!
333   if(!test){
334      for(unsigned int i=0;i<m_vRootFileNames.size();i++){
335         const string rmfiles("rm -f "+m_vRootFileNames[i]);
336         const int res = system(rmfiles.c_str());
337         if(res) {
338                cout<<"Could not remove "<<m_vRootFileNames[i]<<endl;
339                cout<<"Please remove it manually!"<<endl;
340         }
341      }
342      const string rmfiles("rm -f "+dir+"/*");
343      system(rmfiles.c_str());
344      const string rmdir="rm -f -r "+dir;
345      const int res2 = system(rmdir.c_str());
346      if(res2) {
347             cout<<"Could not remove "<<dir<<endl;
348             cout<<"Please remove it manually!"<<endl;
349      }
350   }
351 }
352
353 /*******************************************************************************************/
354 //find out what kind of tree we have and call the right merger
355 bool GateMergeManager::MergeTree(string chainName){
356 if (m_fastMerge==false)
357  {
358    TChain* chain = new TChain(chainName.c_str());
359
360    // number of files to merge
361    int nfiles=m_vRootFileNames.size();
362
363    for(int i=0;i<nfiles;i++) chain->Add(m_vRootFileNames[i].c_str());
364    int nentries=chain->GetEntries();
365    if(nentries<=0) {
366       if(m_verboseLevel>1) cout<<chain->GetName()<<" is empty!"<<endl;
367       return false;
368    }
369
370    if(chainName=="Gate") MergeGate(chain);
371
372    if(chain->FindBranch("eventID1")!=NULL) {
373      if(  (chain->FindBranch("runID")==NULL)
374         ||(chain->FindBranch("time1")==NULL) 
375         ||(chain->FindBranch("time2")==NULL) ) {
376         cout<<"Cannot find one of:  runID, time1, time2 in "<<chain->GetName()<<endl;
377         return false;
378      }
379      // coincidence
380      return MergeCoin(chain);
381
382    } else if(  (chain->FindBranch("runID")==NULL)||(chain->FindBranch("eventID")==NULL)
383                                                   ||(chain->FindBranch("time")==NULL) )
384             {
385              cout<<"Cannot find one of: runID, eventID, time in "<<chain->GetName()<<endl;
386              return false;
387             }
388            // Singles or Hits
389            return MergeSing(chain);
390    }
391    else
392    {
393      if (chainName.find("Hits",0)!=string::npos) return FastMergeSing(chainName);
394      if (chainName.find("Singles",0)!=string::npos) return FastMergeSing(chainName);
395      if (chainName.find("Gate",0)!=string::npos) return FastMergeGate(chainName);
396      if (chainName.find("Coincidences",0)!=string::npos) return FastMergeCoin(chainName);
397    }
398    return 0;
399 }
400
401 /*******************************************************************************************/
402 bool GateMergeManager::FastMergeGate(string name)
403 {
404 //create the output file
405 float   event   = 0;
406 int offset      = 0;
407 int currentfile = 0;
408 float lastEvent =-1;
409 //float maxtime   =-999999999;
410 string clusterName=name+"_cluster";
411
412 for(int j=0;j<m_Nfiles;j++)
413  {
414    filearr[j]->ReOpen("UPDATE");
415    TTree *oldTree = (TTree*)gDirectory->Get(name.c_str());
416    TBranch *branch  = oldTree->GetBranch("event");
417    branch->SetAddress(&event);
418    TBranch* newBranch=oldTree->Branch("eventcluster",&event,"event");
419
420    Int_t nentries=(Int_t)oldTree->GetEntries();
421
422    if (j==0) offset=0;
423    else offset+=m_lastEvents[currentfile];
424    currentfile++;
425
426    for(int i=0;i<nentries;i++)
427       {
428        branch->GetEvent(i);
429        if(lastEvent!=event)
430         {
431          lastEvent=event;
432          event+=offset;
433          newBranch->Fill();
434         } 
435        else
436         {
437          if((i==nentries-1) && (j==m_Nfiles-1))
438           {
439            event+=offset;
440            newBranch->Fill();
441           }
442          else if (i!=nentries-1) 
443               { 
444                branch->GetEvent(i);
445                event+=offset;
446                newBranch->Fill();
447                lastEvent=0;
448                offset=0;
449               }
450         }
451       }
452     oldTree->Write();
453  }
454  return true;
455 }
456
457 /*******************************************************************************************/
458 bool GateMergeManager::FastMergeSing(string name)
459 {
460 int eventID = 0;
461 int runID = 0; 
462 int offset      = 0;
463 int currentfile = 0;
464 float lastRun =-1;
465 string clusterName=name+"_cluster"; 
466
467 for(int j=0;j<m_Nfiles;j++)
468  {
469    filearr[j]->ReOpen("UPDATE");
470    TTree *oldTree = (TTree*)gDirectory->Get(name.c_str());
471    TBranch *branch  = oldTree->GetBranch("eventID");
472    branch->SetAddress(&eventID);
473    TBranch *branch2  = oldTree->GetBranch("runID");
474    branch2->SetAddress(&runID);
475    TBranch* newBranch=oldTree->Branch("eventIDcluster",&eventID,"eventID/I");
476
477    Int_t nentries=(Int_t)oldTree->GetEntries();
478
479    if (j==0) offset=0;
480    else offset+=m_lastEvents[currentfile];
481    currentfile++;
482
483    for(int i=0;i<nentries;i++)
484      {
485       branch->GetEvent(i);
486       if(lastRun!=runID)
487         {
488          lastRun=runID;
489          offset=0;
490         }
491       eventID+=offset;
492       newBranch->Fill();
493     }
494   oldTree->Write();
495   }
496   return true;
497 }
498
499 /*******************************************************************************************/
500 bool GateMergeManager::FastMergeCoin(string name)
501 {
502 //create the output file
503 int eventID1 = 0;
504 int eventID2 = 0;
505 int runID = 0;
506 int offset      = 0;
507 int currentfile = 0;
508 float lastRun =-1;
509 string clusterName=name+"_cluster";
510 cout<<"starting coincidence merging..."<<endl;
511 for(int j=0;j<m_Nfiles;j++)
512  {
513 cout<<"working on file..."<<j<<endl;
514    filearr[j]->ReOpen("UPDATE");
515    TTree *oldTree = (TTree*)gDirectory->Get(name.c_str());
516    TBranch *branch1  = oldTree->GetBranch("eventID1");
517    branch1->SetAddress(&eventID1);
518    TBranch *branch2  = oldTree->GetBranch("eventID2");
519    branch2->SetAddress(&eventID2);
520    TBranch *branch3  = oldTree->GetBranch("runID");
521    branch3->SetAddress(&runID);
522
523    TBranch* newBranch1=oldTree->Branch("eventID1cluster",&eventID1,"eventID1/I");
524    TBranch* newBranch2=oldTree->Branch("eventID2cluster",&eventID2,"eventID2/I");
525
526    Int_t nentries=(Int_t)oldTree->GetEntries();
527
528    if (j==0) offset=0;
529    else offset+=m_lastEvents[currentfile];
530    cout<<"the offset in this file is: "<<offset<<endl;
531    currentfile++;
532
533    for(int i=0;i<nentries;i++){
534        branch1->GetEvent(i);
535        branch2->GetEvent(i);
536        branch3->GetEvent(i);
537        if(lastRun!=runID) 
538         {
539          lastRun=runID;
540          offset=0; 
541         }
542        eventID1+=offset;
543        eventID2+=offset;
544        newBranch1->Fill();
545        newBranch2->Fill();
546     }
547     oldTree->Write();
548  }
549     return true;
550 }
551
552 /*******************************************************************************************/
553 // Gate tree merger
554 bool GateMergeManager::MergeGate(TChain* chainG) {
555
556    int nentries=chainG->GetEntries();   
557
558    float   event   = 0;
559    chainG->SetBranchAddress("event",&event);
560    Float_t iontime = 0;
561    chainG->SetBranchAddress("iontime",&iontime);
562
563    //create the new tree
564    //Allow for large files and do not write every bit separately
565    m_RootTarget->cd();
566    TTree * newTree = chainG->CloneTree(0);
567    newTree->SetAutoSave(2000000000);
568    if(m_maxRoot!=0) newTree->SetMaxTreeSize(m_maxRoot);
569    else newTree->SetMaxTreeSize(17179869184LL);
570
571    // changing the compression level everywhere
572    TBranch *br;
573    TIter next(newTree->GetListOfBranches());
574    while ((br=(TBranch*)next())) br->SetCompressionLevel(m_CompLevel);
575
576    // the main part: modify the eventID
577    int offset      = 0;
578    int currentTree = 0;
579    float lastEvent =-1;
580    float maxtime   =-999999999;
581    for(int i=0;i<nentries;i++){
582        if (chainG->GetEntry(i) <= 0) break; //end of chain
583        if(chainG->GetTreeNumber()>currentTree) {
584            currentTree++;
585            offset+=m_lastEvents[currentTree];
586            // check for overlaping time intervalls between different files
587            // (not within the same file i.e. no time order assumed)
588            if(iontime<maxtime)
589               if(m_verboseLevel>0) 
590                  cout<<"Warning - overlapping Gate iontime ("
591                      <<iontime<<") in file: "<<m_vRootFileNames[currentTree].c_str()<<endl;
592        }
593        // run end is marked by repeating event
594        if(lastEvent!=event) {
595            lastEvent=event;
596 //           if (fpclassify(iontime)!=FP_NAN && fpclassify(iontime)!=FP_INFINITE) {
597            if(finite(iontime)){
598              if(maxtime<iontime) maxtime=iontime;
599              // the offset to get a unique event numbering
600              event+=offset;
601              newTree->Fill();
602            } else {
603              cout<<"Warning - inf or NaN in Gate tree for iontime! "<<m_vRootFileNames[currentTree].c_str()<<endl;
604              //we are lost anyhow
605              maxtime   =-999999999;
606            }
607         } else {
608            // run end 
609            if(chainG->GetEntry(i+1) <= 0)  {              //chain end fill the double event 
610               event+=offset;
611               newTree->Fill();
612            } else if(chainG->GetTreeNumber()==currentTree             //no new file or
613                     ||(chainG->GetTreeNumber()>currentTree
614                        &&m_lastEvents[chainG->GetTreeNumber()]==0) ) { //new file with no offset fill double event
615               chainG->GetEntry(i);
616               event+=offset;
617               newTree->Fill();
618               lastEvent=0;// we may have triple 1 (1 event run)
619               offset=0;
620            }
621         }
622       }
623       newTree->Write();
624       delete newTree;
625       return true;
626 }
627
628 /*******************************************************************************************/
629 // Singles and Hits tree merger
630 bool GateMergeManager::MergeSing(TChain* chainS){
631
632    int nentriesS=chainS->GetEntries();
633
634    int eventID = 0;
635    int runID   = 0;
636    double time = 0;
637
638    chainS->SetBranchAddress("eventID",&eventID);
639    chainS->SetBranchAddress("runID",&runID);
640    chainS->SetBranchAddress("time",&time);
641
642    m_RootTarget->cd();
643    TTree * newSing = chainS->CloneTree(0);
644    newSing->SetAutoSave(2000000000);
645    if(m_maxRoot!=0) newSing->SetMaxTreeSize(m_maxRoot);
646    else newSing->SetMaxTreeSize(17179869184LL);
647
648    // changing CompLevel everywhere
649    TBranch *br;
650    TIter next(newSing->GetListOfBranches());
651    while ((br=(TBranch*)next())) br->SetCompressionLevel(m_CompLevel);
652
653    //the main part: changing eventID
654    int offset      = 0;
655    int currentTree = 0;
656    float lastRun   = 0;
657    float maxtime   =-999999999;
658    for(int i=0;i<nentriesS;i++){
659          if(chainS->GetEntry(i) <= 0) break; //end of chain
660          if(chainS->GetTreeNumber()>currentTree) {
661            currentTree++;
662             offset+=m_lastEvents[currentTree];
663            // check for overlaping time intervalls between different files
664            // (not within the same file i.e. no time order assumed)
665            if(time<maxtime) 
666                if(m_verboseLevel>0) 
667                   cout<<"Warning - overlapping Singles time ("
668                       <<time<<") in file: "<<m_vRootFileNames[currentTree].c_str()<<endl;
669          }
670          // the offset to get a unique event numbering
671          if(lastRun!=runID) {
672             // run end 
673             lastRun=runID;
674             offset=0;     //new run in file we must not change the eventID anymore
675          }
676          eventID+=offset;
677          newSing->Fill();
678          if(maxtime<time)maxtime=time;
679     }
680     newSing->Write();
681     delete newSing;
682     return true;
683 };
684
685 /*******************************************************************************************/
686 // Coincidences tree merger
687 bool GateMergeManager::MergeCoin(TChain* chainC){
688
689    int nentriesC=chainC->GetEntries();
690
691     int eventID1 = 0;
692     chainC->SetBranchAddress("eventID1",&eventID1);
693     double time1 = 0;
694     chainC->SetBranchAddress("time1",&time1);
695     int eventID2 = 0;
696     chainC->SetBranchAddress("eventID2",&eventID2);
697     double time2 = 0;
698     chainC->SetBranchAddress("time2",&time2);
699     int runID    = 0;
700     chainC->SetBranchAddress("runID",&runID);
701
702     m_RootTarget->cd();
703     TTree * newCoin = chainC->CloneTree(0);
704     newCoin->SetAutoSave(2000000000);
705     if(m_maxRoot!=0) newCoin->SetMaxTreeSize(m_maxRoot);
706     else newCoin->SetMaxTreeSize(17179869184LL);
707
708     // changing CompLevel everywhere
709     TBranch *br;
710     TIter next(newCoin->GetListOfBranches());
711     while ((br=(TBranch*)next())) br->SetCompressionLevel(m_CompLevel);
712
713     int offset      = 0;
714     int currentTree = 0;
715     float lastRun   = 0;
716     float maxtime   =-999999999;
717     for(int i=0;i<nentriesC;i++){
718        if(chainC->GetEntry(i) <= 0) break; //end of chain
719        if(chainC->GetTreeNumber()>currentTree) {
720          currentTree++;
721          offset+=m_lastEvents[currentTree];
722          // check for overlaping time intervalls between different files
723          // (not within the same file i.e. no time order assumed)
724          if(time1<maxtime) 
725            if(m_verboseLevel>0)
726                cout<<"Warning - overlapping Coincidences time1 ("
727                    <<time1<<") in file: "<<m_vRootFileNames[currentTree].c_str()<<endl;
728          if(time2<maxtime) 
729            if(m_verboseLevel>0)
730                cout<<"Warning - overlapping Coincidences time2 ("
731                    <<time2<<") in file: "<<m_vRootFileNames[currentTree].c_str()<<endl;
732        }
733        // the offset to get a unique event numbering
734        if(lastRun!=runID) {
735           // run end
736           lastRun=runID;
737           offset=0; //new run in file we must not change the eventID anymore
738        }
739        eventID1+=offset;
740        eventID2+=offset;
741        newCoin->Fill();
742        if(maxtime<time1)maxtime=time1;
743        if(maxtime<time2)maxtime=time2;
744     }
745     newCoin->Write();
746     delete newCoin;
747     return true;
748 }
749 /*******************************************************************************************/