TTreeProxyGenerator.cxx

Go to the documentation of this file.
00001 // @(#)root/treeplayer:$Id: TTreeProxyGenerator.cxx 36449 2010-10-28 20:52:17Z pcanal $
00002 // Author: Philippe Canal 06/06/2004
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers and al.        *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 /*
00013   TODO:
00014   Have separate names for the wrapper classes in the cases of: [done]
00015   clones/non clones
00016   split/non split
00017   split levels
00018 
00019   Have a solution for passing top+"."+middle to the parents classes [probably done .. need testing]
00020 
00021   Have a solution for the return by references of abstract classes [not done]
00022 
00023   Have object inside ClonesArray properly treated! [done]
00024   Why is there 2 TRef proxy classes? [done]
00025 
00026   check why some inheritance are TObjProxy and not TPx_
00027 
00028   Be smart enough to avoid issue about having 2 classes one unrolled and one non unrolled!
00029 
00030   When using in interpreted mode understand why the reloading reloads the calling script and then crashes :(
00031 
00032   CINT does not properly call the custom operators when doing return fNtrack.
00033 
00034   CINT does not handle fMatrix[2][1] well.
00035 
00036   The user's function in script.h are not exposed by ACLiC.
00037 
00038   Review the method to avoid the useless refreshing of the generated file
00039   - for most efficiency it would require a different name for each tree
00040 */
00041 
00042 #include "TTreeProxyGenerator.h"
00043 
00044 #include "TFriendProxyDescriptor.h"
00045 #include "TBranchProxyDescriptor.h"
00046 #include "TBranchProxyClassDescriptor.h"
00047 
00048 #include "TList.h"
00049 #include "Varargs.h"
00050 #include <stdio.h>
00051 
00052 class TTree;
00053 class TBranch;
00054 class TStreamerElement;
00055 
00056 #include "TClass.h"
00057 #include "TClassEdit.h"
00058 #include "TClonesArray.h"
00059 #include "TError.h"
00060 #include "TROOT.h"
00061 #include "TObjString.h"
00062 
00063 #include "TTreeFormula.h"
00064 #include "TFormLeafInfo.h"
00065 
00066 
00067 #include "TBranchElement.h"
00068 #include "TChain.h"
00069 #include "TFile.h"
00070 #include "TFriendElement.h"
00071 #include "TLeaf.h"
00072 #include "TTree.h"
00073 #include "TVirtualStreamerInfo.h"
00074 #include "TStreamerElement.h"
00075 #include "TSystem.h"
00076 #include "TLeafObject.h"
00077 #include "TVirtualCollectionProxy.h"
00078 
00079 void Debug(Int_t level, const char *va_(fmt), ...)
00080 {
00081    // Use this function in case an error occured.
00082 
00083    if (gDebug>=level) {
00084       va_list ap;
00085       va_start(ap,va_(fmt));
00086       ErrorHandler(kInfo,"TTreeProxyGenerator",va_(fmt), ap);
00087       va_end(ap);
00088    }
00089 }
00090 
00091 namespace {
00092 
00093    Bool_t AreDifferent(const TString& from, const TString& to)
00094    {
00095       FILE *left = fopen(from.Data(),"r");
00096       FILE *right = fopen(to.Data(),"r");
00097 
00098       char leftbuffer[256];
00099       char rightbuffer[256];
00100 
00101       char *lvalue,*rvalue;
00102 
00103       Bool_t areEqual = kTRUE;
00104 
00105       do {
00106          lvalue = fgets(leftbuffer, sizeof(leftbuffer), left);
00107          rvalue = fgets(rightbuffer, sizeof(rightbuffer), right);
00108 
00109          if (lvalue&&rvalue) {
00110             if (strstr(lvalue,"by ROOT version")) {
00111                // skip the comment line with the time and date
00112             } else {
00113                areEqual = areEqual && (0 == strncmp(lvalue,rvalue,sizeof(leftbuffer)));
00114             }
00115          }
00116          if (lvalue&&!rvalue) areEqual = kFALSE;
00117          if (rvalue&&!lvalue) areEqual = kFALSE;
00118 
00119       } while(areEqual && lvalue && rvalue);
00120 
00121       fclose(left);
00122       fclose(right);
00123 
00124       return !areEqual;
00125    }
00126 }
00127 
00128 namespace ROOT {
00129 
00130    TString GetArrayType(TStreamerElement *element, const char *subtype,
00131                         TTreeProxyGenerator::EContainer container)
00132    {
00133       TString result;
00134       int ndim = 0;
00135       if (element->InheritsFrom(TStreamerBasicPointer::Class())) {
00136          TStreamerBasicPointer * elem = (TStreamerBasicPointer*)element;
00137          const char *countname = elem->GetCountName();
00138          if (countname && strlen(countname)>0) ndim = 1;
00139       }
00140       ndim += element->GetArrayDim();
00141 
00142       TString middle;
00143       if (container == TTreeProxyGenerator::kClones) {
00144          middle = "Cla";
00145       } else if  (container == TTreeProxyGenerator::kSTL) {
00146          middle = "Stl";
00147       }
00148 
00149       if (ndim==0) {
00150          result = "T";
00151          result += middle;
00152          result += subtype;
00153          result += "Proxy";
00154       } else if (ndim==1) {
00155          result = "T";
00156          result += middle;
00157          result += "Array";
00158          result += subtype;
00159          result += "Proxy";
00160       } else {
00161          result = "T";
00162          result += middle;
00163          result += "ArrayProxy<";
00164          for(Int_t ind = ndim - 2; ind > 0; --ind) {
00165             result += "TMultiArrayType<";
00166          }
00167          result += "TArrayType<";
00168          result += element->GetTypeName();
00169          result += ",";
00170          result += element->GetMaxIndex(ndim-1);
00171          result += "> ";
00172          for(Int_t ind = ndim - 2; ind > 0; --ind) {
00173             result += ",";
00174             result += element->GetMaxIndex(ind);
00175             result += "> ";
00176          }
00177          result += ">";
00178       }
00179       return result;
00180 
00181       /*
00182         if (!strcmp("unsigned int", name))
00183         sprintf(line, "%u", *(unsigned int *)buf);
00184         else if (!strcmp("int", name))
00185         sprintf(line, "%d", *(int *)buf);
00186         else if (!strcmp("unsigned long", name))
00187         sprintf(line, "%lu", *(unsigned long *)buf);
00188         else if (!strcmp("long", name))
00189         sprintf(line, "%ld", *(long *)buf);
00190         else if (!strcmp("unsigned short", name))
00191         sprintf(line, "%hu", *(unsigned short *)buf);
00192         else if (!strcmp("short", name))
00193         sprintf(line, "%hd", *(short *)buf);
00194         else if (!strcmp("unsigned char", name))
00195         sprintf(line, "%u", *(unsigned char *)buf);
00196         else if (!strcmp("bool", name))
00197         sprintf(line, "%u", *(unsigned char *)buf);
00198         else if (!strcmp("char", name))
00199         sprintf(line, "%d", *(char *)buf);
00200         else if (!strcmp("float", name))
00201         sprintf(line, "%g", *(float *)buf);
00202         else if (!strcmp("double", name))
00203         sprintf(line, "%g", *(double *)buf);
00204       */
00205    }
00206 
00207    TTreeProxyGenerator::TTreeProxyGenerator(TTree* tree,
00208                                             const char *script,
00209                                             const char *fileprefix,
00210                                             const char *option, UInt_t maxUnrolling) :
00211       fMaxDatamemberType(2),
00212       fScript(script),
00213       fCutScript(),
00214       fPrefix(fileprefix),
00215       fHeaderFileName(),
00216       fOptionStr(option),
00217       fOptions(0),
00218       fMaxUnrolling(maxUnrolling),
00219       fTree(tree),
00220       fCurrentListOfTopProxies(&fListOfTopProxies)
00221    {
00222       // Constructor.
00223 
00224       ParseOptions();
00225 
00226       AnalyzeTree(fTree);
00227 
00228       WriteProxy();
00229    }
00230 
00231    TTreeProxyGenerator::TTreeProxyGenerator(TTree* tree,
00232                                             const char *script, const char *cutscript,
00233                                             const char *fileprefix,
00234                                             const char *option, UInt_t maxUnrolling) :
00235       fMaxDatamemberType(2),
00236       fScript(script),
00237       fCutScript(cutscript),
00238       fPrefix(fileprefix),
00239       fHeaderFileName(),
00240       fOptionStr(option),
00241       fOptions(0),
00242       fMaxUnrolling(maxUnrolling),
00243       fTree(tree),
00244       fCurrentListOfTopProxies(&fListOfTopProxies)
00245    {
00246       // Constructo.
00247 
00248       ParseOptions();
00249 
00250       AnalyzeTree(fTree);
00251 
00252       WriteProxy();
00253    }
00254 
00255    Bool_t TTreeProxyGenerator::NeedToEmulate(TClass *cl, UInt_t /* level */)
00256    {
00257       // Return true if we should create a nested class representing this class
00258       
00259       return cl!=0 && cl->TestBit(TClass::kIsEmulation);
00260    }
00261 
00262    TBranchProxyClassDescriptor*
00263    TTreeProxyGenerator::AddClass( TBranchProxyClassDescriptor* desc )
00264    {
00265       // Add a Class Descriptor.
00266 
00267       if (desc==0) return 0;
00268 
00269       TBranchProxyClassDescriptor *existing =
00270          (TBranchProxyClassDescriptor*)fListOfClasses(desc->GetName());
00271 
00272       int count = 0;
00273       while (existing) {
00274          if (! existing->IsEquivalent( desc )  ) {
00275             TString newname = desc->GetRawSymbol();
00276             count++;
00277             newname += "_";
00278             newname += count;
00279 
00280             desc->SetName(newname);
00281             existing = (TBranchProxyClassDescriptor*)fListOfClasses(desc->GetName());
00282          } else {
00283             // we already have the exact same class
00284             delete desc;
00285             return existing;
00286          }
00287       }
00288       fListOfClasses.Add(desc);
00289       return desc;
00290    }
00291 
00292    void TTreeProxyGenerator::AddFriend( TFriendProxyDescriptor* desc )
00293    {
00294       // Add Friend descriptor.
00295 
00296       if (desc==0) return;
00297 
00298       TFriendProxyDescriptor *existing =
00299          (TFriendProxyDescriptor*)fListOfFriends(desc->GetName());
00300 
00301       int count = 0;
00302       while (existing) {
00303          if (! existing->IsEquivalent( desc )  ) {
00304             TString newname = desc->GetName();
00305             count++;
00306             newname += "_";
00307             newname += count;
00308 
00309             desc->SetName(newname);
00310             existing = (TFriendProxyDescriptor*)fListOfFriends(desc->GetName());
00311 
00312          } else {
00313 
00314             desc->SetDuplicate();
00315             break;
00316          }
00317       }
00318 
00319       // Insure uniqueness of the title also.
00320       TString basetitle = desc->GetTitle();
00321       TIter next( &fListOfFriends );
00322       while ( (existing = (TFriendProxyDescriptor*)next()) ) {
00323          if (strcmp(existing->GetTitle(),desc->GetTitle())==0) {
00324 
00325             TString newtitle = basetitle;
00326             count++;
00327             newtitle += "_";
00328             newtitle += count;
00329 
00330             desc->SetTitle(newtitle);
00331 
00332             // Restart of the begining of the loop.
00333             next = &fListOfFriends;
00334          }
00335       }
00336 
00337       fListOfFriends.Add(desc);
00338    }
00339 
00340    void TTreeProxyGenerator::AddForward( const char *classname )
00341    {
00342       // Add a forward declaration request.
00343 
00344       TObject *obj = fListOfForwards.FindObject(classname);
00345       if (obj) return;
00346 
00347       if (strstr(classname,"<")!=0) {
00348          // this is a template instantiation.
00349          // let's ignore it for now
00350 
00351          if (gDebug>=6) Warning("AddForward","Forward declaration of templated class not implemented yet.");
00352       } else if (strcmp(classname,"string")==0) {
00353          // no need to forward declare string
00354       } else {
00355          fListOfForwards.Add(new TNamed(classname,Form("class %s;\n",classname)));
00356       }
00357       return;
00358    }
00359 
00360    void TTreeProxyGenerator::AddForward(TClass *cl)
00361    {
00362       // Add a forward declaration request.
00363 
00364       if (cl) AddForward(cl->GetName());
00365    }
00366 
00367    void TTreeProxyGenerator::AddHeader(TClass *cl)
00368    {
00369       // Add a header inclusion request.
00370 
00371       if (cl==0) return;
00372 
00373       TObject *obj = fListOfHeaders.FindObject(cl->GetName());
00374       if (obj) return;
00375 
00376       TString directive;
00377 
00378       if (cl->GetCollectionProxy() && cl->GetCollectionProxy()->GetValueClass()) {
00379          AddHeader( cl->GetCollectionProxy()->GetValueClass() );
00380       }
00381       Int_t stlType;
00382       if (cl->GetCollectionProxy() && (stlType=TClassEdit::IsSTLCont(cl->GetName()))) {
00383          const char *what = "";
00384          switch(stlType)  {
00385             case TClassEdit::kVector:   what = "vector"; break;
00386             case TClassEdit::kList:     what = "list"; break;
00387             case -TClassEdit::kDeque: // same as positive
00388             case TClassEdit::kDeque:    what = "deque"; break;
00389             case -TClassEdit::kMap: // same as positive
00390             case TClassEdit::kMap:      what = "map"; break;
00391             case -TClassEdit::kMultiMap: // same as positive
00392             case TClassEdit::kMultiMap: what = "map"; break;
00393             case -TClassEdit::kSet:  // same as positive
00394             case TClassEdit::kSet:      what = "set"; break;
00395             case -TClassEdit::kMultiSet: // same as positive
00396             case TClassEdit::kMultiSet: what = "set"; break;
00397          }
00398          if (what[0]) {
00399             directive = "#include <";
00400             directive.Append(what);
00401             directive.Append(">\n");
00402          }
00403       } else if (cl->GetDeclFileName() && strlen(cl->GetDeclFileName()) ) {
00404          // Actually we probably should look for the file ..
00405          const char *filename = cl->GetDeclFileName();
00406 
00407          if (!filename) return;
00408 
00409 #ifdef R__WIN32
00410          TString inclPath("include;prec_stl"); // GetHtml()->GetIncludePath());
00411 #else
00412          TString inclPath("include:prec_stl"); // GetHtml()->GetIncludePath());
00413 #endif
00414          Ssiz_t posDelim = 0;
00415          TString inclDir;
00416          TString sIncl(filename);
00417 #ifdef R__WIN32
00418          const char* pdelim = ";";
00419          static const char ddelim = '\\';
00420 #else
00421          const char* pdelim = ":";
00422          static const char ddelim = '/';
00423 #endif
00424          while (inclPath.Tokenize(inclDir, posDelim, pdelim))
00425          {
00426             if (sIncl.BeginsWith(inclDir)) {
00427                filename += inclDir.Length();
00428                if (filename[0] == ddelim || filename[0] == '/') {
00429                   ++filename;
00430                }
00431                break;
00432             }
00433          }
00434          directive = Form("#include \"%s\"\n",filename);
00435       } else if (!strncmp(cl->GetName(), "pair<", 5)
00436                  || !strncmp(cl->GetName(), "std::pair<", 10)) {
00437          TClassEdit::TSplitType split(cl->GetName());
00438          if (split.fElements.size() == 3) {
00439             for (int arg = 1; arg < 3; ++arg) {
00440                TClass* clArg = TClass::GetClass(split.fElements[arg].c_str());
00441                if (clArg) AddHeader(clArg);
00442             }
00443          }
00444       }
00445       if (directive.Length()) {
00446          TIter i( &fListOfHeaders );
00447          for(TNamed *n = (TNamed*) i(); n; n = (TNamed*)i() ) {
00448             if (directive == n->GetTitle()) {
00449                return;
00450             }
00451          }
00452          fListOfHeaders.Add(new TNamed(cl->GetName(),directive.Data()));
00453       }
00454    }
00455 
00456    void TTreeProxyGenerator::AddHeader(const char *classname)
00457    {
00458       // Add a header inclusion request.
00459 
00460       AddHeader(TClass::GetClass(classname));
00461    }
00462 
00463    void TTreeProxyGenerator::AddPragma(const char *pragma_text)
00464    {
00465       // Add a forward declaration request.
00466 
00467       TIter i( &fListOfPragmas );
00468       for(TObjString *n = (TObjString*) i(); n; n = (TObjString*)i() ) {
00469          if (pragma_text == n->GetString()) {
00470             return;
00471          }
00472       }
00473 
00474       fListOfPragmas.Add( new TObjString( pragma_text ) );
00475 
00476    }
00477 
00478    void TTreeProxyGenerator::AddDescriptor(TBranchProxyDescriptor *desc)
00479    {
00480       // Add a branch descriptor.
00481 
00482       if (desc) {
00483          TBranchProxyDescriptor *existing =
00484             (TBranchProxyDescriptor*)((*fCurrentListOfTopProxies)(desc->GetName()));
00485          if (existing) {
00486             Warning("TTreeProxyGenerator","The branch name \"%s\" is duplicated. Only the first instance \n"
00487                "\twill be available directly. The other instance(s) might be available via their complete name\n"
00488                "\t(including the name of their mother branche's name).",desc->GetName());
00489          } else {
00490             fCurrentListOfTopProxies->Add(desc);
00491             UInt_t len = strlen(desc->GetTypeName());
00492             if ((len+2)>fMaxDatamemberType) fMaxDatamemberType = len+2;
00493          }
00494       }
00495    }
00496 
00497 static TString GetContainedClassName(TBranchElement *branch, TStreamerElement *element, Bool_t ispointer)
00498 {
00499    TString cname = branch->GetClonesName();
00500    if (cname.Length()==0) {
00501       // We may have any unsplit clones array
00502       Long64_t i = branch->GetTree()->GetReadEntry();
00503       if (i<0) i = 0;
00504       branch->GetEntry(i);
00505       char *obj = branch->GetObject();
00506 
00507 
00508       TBranchElement *parent = (TBranchElement*)branch->GetMother()->GetSubBranch(branch);
00509       const char *pclname = parent->GetClassName();
00510 
00511       TClass *clparent = TClass::GetClass(pclname);
00512       // TClass *clm = TClass::GetClass(GetClassName());
00513       Int_t lOffset = 0; // offset in the local streamerInfo.
00514       if (clparent) {
00515          const char *ename = 0;
00516          if (element) {
00517             ename = element->GetName();
00518             lOffset = clparent->GetStreamerInfo()->GetOffset(ename);
00519          } else {
00520             lOffset = 0;
00521          }
00522       }
00523       else Error("AnalyzeBranch", "Missing parent for %s.", branch->GetName());
00524 
00525       TClonesArray *arr;
00526       if (ispointer) {
00527          arr = (TClonesArray*)*(void**)(obj+lOffset);
00528       } else {
00529          arr = (TClonesArray*)(obj+lOffset);
00530       }
00531       cname = arr->GetClass()->GetName();
00532 
00533    }
00534    if (cname.Length()==0) {
00535       Error("AnalyzeBranch",
00536          "Introspection of TClonesArray in older file not implemented yet.");
00537    }
00538    return cname;
00539 }
00540 
00541 static TVirtualStreamerInfo *GetStreamerInfo(TBranch *branch, TIter current, TClass *cl)
00542 {
00543    // Return the correct TStreamerInfo of class 'cname' in the list of
00544    // branch (current) [Assuming these branches correspond to a flattened
00545    // version of the class.]
00546 
00547    TVirtualStreamerInfo *objInfo = 0;
00548    TBranchElement *b = 0;
00549    TString cname = cl->GetName();
00550 
00551    while( ( b = (TBranchElement*)current() ) ) {
00552       if ( cname == b->GetInfo()->GetName() ) {
00553          objInfo = b->GetInfo();
00554          break;
00555       }
00556    }
00557    if (objInfo == 0 && branch->GetTree()->GetDirectory()->GetFile()) {
00558       TVirtualStreamerInfo *i = (TVirtualStreamerInfo *)branch->GetTree()->GetDirectory()->GetFile()->GetStreamerInfoCache()->FindObject(cname);
00559       if (i) {
00560          // NOTE: Is this correct for Foreigh classes?
00561          objInfo = (TVirtualStreamerInfo *)cl->GetStreamerInfo(i->GetClassVersion());
00562       }
00563    }
00564    if (objInfo == 0) {
00565       // We still haven't found it ... this is likely to be an STL collection .. anyway, use the current StreamerInfo.
00566       objInfo = cl->GetStreamerInfo();
00567    }
00568    return objInfo;
00569 }
00570 
00571 static TVirtualStreamerInfo *GetBaseClass(TStreamerElement *element)
00572 {
00573    TStreamerBase *base = dynamic_cast<TStreamerBase*>(element);
00574    if (base) {
00575       TVirtualStreamerInfo *info = base->GetClassPointer()->GetStreamerInfo(base->GetBaseVersion());
00576       if (info) return info;
00577    }
00578    return 0;
00579 }
00580 
00581    UInt_t TTreeProxyGenerator::AnalyzeBranches(UInt_t level,TBranchProxyClassDescriptor *topdesc,
00582                                                TBranchElement *branch, TVirtualStreamerInfo *info)
00583    {
00584       // Analyze the sub-branch and populate the TTreeProxyGenerator or the topdesc with
00585       // its findings.
00586 
00587       if (info==0) info = branch->GetInfo();
00588 
00589       TIter branches( branch->GetListOfBranches() );
00590 
00591       return AnalyzeBranches( level, topdesc, branches, info );
00592    }
00593 
00594    UInt_t TTreeProxyGenerator::AnalyzeBranches(UInt_t level,
00595                                                TBranchProxyClassDescriptor *topdesc,
00596                                                TIter &branches,
00597                                                TVirtualStreamerInfo *info)
00598    {
00599       // Analyze the list of sub branches of a TBranchElement by looping over
00600       // the streamer elements and create the appropriate class proxies.
00601 
00602 /*
00603 
00604    Find the content class name (GetClassName)
00605    Record wether this is a collection or not
00606 
00607    Find the StreamerInfo
00608 
00609    For each streamerelement
00610       if element is base
00611          if name match, loop over subbranches
00612          otherwise loop over current branches
00613       else if eleement is object (or pointer to object?)
00614          if name match go ahead, loop over subbranches
00615          if name does not match. loop over current branches (fix names).
00616       else
00617          add branch.
00618 
00619 */
00620       UInt_t lookedAt = 0;
00621       EContainer container = kNone;
00622       TString middle;
00623       TString proxyTypeName;
00624       TBranchProxyClassDescriptor::ELocation outer_isclones = TBranchProxyClassDescriptor::kOut;
00625       TString containerName;
00626       TString subBranchPrefix;
00627       Bool_t skipped = false;
00628 
00629       {
00630          TIter peek = branches;
00631          TBranchElement *branch = (TBranchElement*)peek();
00632          if (topdesc && topdesc->IsClones()) {
00633             container = kClones;
00634             middle = "Cla";
00635             outer_isclones = TBranchProxyClassDescriptor::kClones;
00636             containerName = "TClonesArray";
00637          } else if (topdesc && topdesc->IsSTL()) {
00638             container = kSTL;
00639             middle = "Stl";
00640             outer_isclones = TBranchProxyClassDescriptor::kSTL;
00641             containerName = topdesc->GetContainerName();
00642          } else if (!topdesc && branch && branch->GetBranchCount() == branch->GetMother()) {
00643             if ( ((TBranchElement*)(branch->GetMother()))->GetType()==3)  {
00644                container = kClones;
00645                middle = "Cla";
00646                outer_isclones = TBranchProxyClassDescriptor::kClones;
00647                containerName = "TClonesArray";
00648             } else {
00649                container = kSTL;
00650                middle = "Stl";
00651                outer_isclones = TBranchProxyClassDescriptor::kSTL;
00652                containerName = branch->GetMother()->GetClassName();
00653             }
00654          } else if (branch->GetType() == 3) {
00655             outer_isclones = TBranchProxyClassDescriptor::kClones;
00656             containerName = "TClonesArray";
00657          } else if (branch->GetType() == 4) {
00658             outer_isclones = TBranchProxyClassDescriptor::kSTL;
00659             containerName = branch->GetMother()->GetSubBranch(branch)->GetClassName();
00660          }
00661          if (topdesc) {
00662             subBranchPrefix = topdesc->GetSubBranchPrefix();
00663          } else {
00664             TBranchElement *mom = (TBranchElement*)branch->GetMother();
00665             subBranchPrefix = mom->GetName();
00666             if (subBranchPrefix[subBranchPrefix.Length()-1]=='.') {
00667                subBranchPrefix.Remove(subBranchPrefix.Length()-1);
00668             } else if (mom->GetType()!=3 && mom->GetType() != 4) {
00669                subBranchPrefix = "";
00670             }
00671          }
00672       }
00673       TIter elements( info->GetElements() );
00674       for( TStreamerElement *element = (TStreamerElement*)elements();
00675            element;
00676            element = (TStreamerElement*)elements() )
00677       {
00678          Bool_t isBase = false;
00679          Bool_t usedBranch = kTRUE;
00680          TString prefix;
00681          TIter peek = branches;
00682          TBranchElement *branch = (TBranchElement*)peek();
00683 
00684          if (branch==0) {
00685             if (topdesc) {
00686                Error("AnalyzeBranches","Ran out of branches when looking in branch %s, class %s",
00687                      topdesc->GetBranchName(), info->GetName());
00688             } else {
00689                Error("AnalyzeBranches","Ran out of branches when looking in class %s, element %s",
00690                      info->GetName(),element->GetName());
00691             }
00692             return lookedAt;
00693          }
00694 
00695          if (info->GetClass()->GetCollectionProxy() && strcmp(element->GetName(),"This")==0) {
00696             // Skip the artifical streamer element.
00697             continue;
00698          }
00699 
00700          if (element->GetType()==-1) {
00701             // This is an ignored TObject base class.
00702             continue;
00703          }
00704 
00705          TString branchname = branch->GetName();
00706          TString branchEndname;
00707          {
00708             TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
00709             if (leaf && outer_isclones == TBranchProxyClassDescriptor::kOut
00710                 && !(branch->GetType() == 3 || branch->GetType() == 4)) branchEndname = leaf->GetName();
00711             else branchEndname = branch->GetName();
00712             Int_t pos;
00713             pos = branchEndname.Index(".");
00714             if (pos!=-1) {
00715                if (subBranchPrefix.Length() &&
00716                   branchEndname.BeginsWith( subBranchPrefix ) ) {
00717                      // brprefix += topdesc->GetSubBranchPrefix();
00718                   branchEndname.Remove(0,subBranchPrefix.Length()+1);
00719                }
00720             }
00721          }
00722 
00723          Bool_t ispointer = false;
00724          switch(element->GetType()) {
00725 
00726             case TVirtualStreamerInfo::kBool:    { proxyTypeName = "T" + middle + "BoolProxy"; break; }
00727             case TVirtualStreamerInfo::kChar:    { proxyTypeName = "T" + middle + "CharProxy"; break; }
00728             case TVirtualStreamerInfo::kShort:   { proxyTypeName = "T" + middle + "ShortProxy"; break; }
00729             case TVirtualStreamerInfo::kInt:     { proxyTypeName = "T" + middle + "IntProxy"; break; }
00730             case TVirtualStreamerInfo::kLong:    { proxyTypeName = "T" + middle + "LongProxy"; break; }
00731             case TVirtualStreamerInfo::kLong64:  { proxyTypeName = "T" + middle + "Long64Proxy"; break; }
00732             case TVirtualStreamerInfo::kFloat:   { proxyTypeName = "T" + middle + "FloatProxy"; break; }
00733             case TVirtualStreamerInfo::kFloat16: { proxyTypeName = "T" + middle + "Float16Proxy"; break; }
00734             case TVirtualStreamerInfo::kDouble:  { proxyTypeName = "T" + middle + "DoubleProxy"; break; }
00735             case TVirtualStreamerInfo::kDouble32:{ proxyTypeName = "T" + middle + "Double32Proxy"; break; }
00736             case TVirtualStreamerInfo::kUChar:   { proxyTypeName = "T" + middle + "UCharProxy"; break; }
00737             case TVirtualStreamerInfo::kUShort:  { proxyTypeName = "T" + middle + "UShortProxy"; break; }
00738             case TVirtualStreamerInfo::kUInt:    { proxyTypeName = "T" + middle + "UIntProxy"; break; }
00739             case TVirtualStreamerInfo::kULong:   { proxyTypeName = "T" + middle + "ULongProxy"; break; }
00740             case TVirtualStreamerInfo::kULong64: { proxyTypeName = "T" + middle + "ULong64Proxy"; break; }
00741             case TVirtualStreamerInfo::kBits:    { proxyTypeName = "T" + middle + "UIntProxy"; break; }
00742 
00743             case TVirtualStreamerInfo::kCharStar: { proxyTypeName = GetArrayType(element,"Char",container); break; }
00744 
00745                // array of basic types  array[8]
00746             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kBool:    { proxyTypeName = GetArrayType(element,"Bool",container ); break; }
00747             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kChar:    { proxyTypeName = GetArrayType(element,"Char",container ); break; }
00748             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kShort:   { proxyTypeName = GetArrayType(element,"Short",container ); break; }
00749             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kInt:     { proxyTypeName = GetArrayType(element,"Int",container ); break; }
00750             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kLong:    { proxyTypeName = GetArrayType(element,"Long",container ); break; }
00751             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kLong64:  { proxyTypeName = GetArrayType(element,"Long64",container ); break; }
00752             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kFloat:   { proxyTypeName = GetArrayType(element,"Float",container ); break; }
00753             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kFloat16: { proxyTypeName = GetArrayType(element,"Float16",container ); break; }
00754             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kDouble:  { proxyTypeName = GetArrayType(element,"Double",container ); break; }
00755             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kDouble32:{ proxyTypeName = GetArrayType(element,"Double32",container ); break; }
00756             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUChar:   { proxyTypeName = GetArrayType(element,"UChar",container ); break; }
00757             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUShort:  { proxyTypeName = GetArrayType(element,"UShort",container ); break; }
00758             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUInt:    { proxyTypeName = GetArrayType(element,"UInt",container ); break; }
00759             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kULong:   { proxyTypeName = GetArrayType(element,"ULong",container ); break; }
00760             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kULong64: { proxyTypeName = GetArrayType(element,"ULong64",container ); break; }
00761             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kBits:    { proxyTypeName = GetArrayType(element,"UInt",container ); break; }
00762 
00763                // pointer to an array of basic types  array[n]
00764             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kBool:    { proxyTypeName = GetArrayType(element,"Bool",container ); break; }
00765             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kChar:    { proxyTypeName = GetArrayType(element,"Char",container ); break; }
00766             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kShort:   { proxyTypeName = GetArrayType(element,"Short",container ); break; }
00767             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kInt:     { proxyTypeName = GetArrayType(element,"Int",container ); break; }
00768             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kLong:    { proxyTypeName = GetArrayType(element,"Long",container ); break; }
00769             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kLong64:  { proxyTypeName = GetArrayType(element,"Long64",container ); break; }
00770             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kFloat:   { proxyTypeName = GetArrayType(element,"Float",container ); break; }
00771             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kFloat16: { proxyTypeName = GetArrayType(element,"Float16",container ); break; }
00772             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kDouble:  { proxyTypeName = GetArrayType(element,"Double",container ); break; }
00773             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kDouble32:{ proxyTypeName = GetArrayType(element,"Double32",container ); break; }
00774             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUChar:   { proxyTypeName = GetArrayType(element,"UChar",container ); break; }
00775             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUShort:  { proxyTypeName = GetArrayType(element,"UShort",container ); break; }
00776             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUInt:    { proxyTypeName = GetArrayType(element,"UInt",container ); break; }
00777             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kULong:   { proxyTypeName = GetArrayType(element,"ULong",container ); break; }
00778             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kULong64: { proxyTypeName = GetArrayType(element,"ULong64",container ); break; }
00779             case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kBits:    { proxyTypeName = GetArrayType(element,"UInt",container ); break; }
00780 
00781                // array counter //[n]
00782             case TVirtualStreamerInfo::kCounter: { proxyTypeName = "T" + middle + "IntProxy"; break; }
00783 
00784 
00785             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObjectp:
00786             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObjectP:
00787             case TVirtualStreamerInfo::kObjectp:
00788             case TVirtualStreamerInfo::kObjectP:
00789             case TVirtualStreamerInfo::kAnyp:
00790             case TVirtualStreamerInfo::kAnyP:
00791             case TVirtualStreamerInfo::kSTL + TVirtualStreamerInfo::kObjectp:
00792             case TVirtualStreamerInfo::kSTL + TVirtualStreamerInfo::kObjectP:
00793             // set as pointers and fall through to the next switches
00794                ispointer = true;
00795             case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObject:
00796             case TVirtualStreamerInfo::kObject:
00797             case TVirtualStreamerInfo::kTString:
00798             case TVirtualStreamerInfo::kTNamed:
00799             case TVirtualStreamerInfo::kTObject:
00800             case TVirtualStreamerInfo::kAny:
00801             case TVirtualStreamerInfo::kBase:
00802             case TVirtualStreamerInfo::kSTL: {
00803                TClass *cl = element->GetClassPointer();
00804                R__ASSERT(cl);
00805 
00806                proxyTypeName = Form("T%sObjProxy<%s >", middle.Data(), cl->GetName());
00807                TString cname = cl->GetName();
00808                TBranchProxyClassDescriptor::ELocation isclones = outer_isclones;
00809                if (cl==TClonesArray::Class()) {
00810                   isclones = TBranchProxyClassDescriptor::kClones;
00811                   cname = GetContainedClassName(branch, element, ispointer);
00812                   containerName = "TClonesArray";
00813                } else if (cl->GetCollectionProxy()) {
00814                   isclones = TBranchProxyClassDescriptor::kSTL;
00815                   containerName = cl->GetName();
00816                   TClass *valueClass = cl->GetCollectionProxy()->GetValueClass();
00817                   if (valueClass) cname = valueClass->GetName();
00818                   else {
00819                      proxyTypeName = Form("TStlSimpleProxy<%s >", cl->GetName());
00820 //                   AddPragma(Form("#pragma create TClass %s;\n", cl->GetName()));
00821                      if (!cl->IsLoaded()) AddPragma(Form("#pragma link C++ class %s;\n", cl->GetName()));
00822                   }
00823                }
00824 
00825                TBranch *parent = branch->GetMother()->GetSubBranch(branch);
00826                TVirtualStreamerInfo *objInfo = 0;
00827                if (branch->GetListOfBranches()->GetEntries()) {
00828                   objInfo = ((TBranchElement*)branch->GetListOfBranches()->At(0))->GetInfo();
00829                } else {
00830                   objInfo = branch->GetInfo();
00831                }
00832                if (element->IsBase()) {
00833                   isBase = true;
00834                   prefix  = "base";
00835 
00836                   if (cl == TObject::Class() && info->GetClass()->CanIgnoreTObjectStreamer())
00837                   {
00838                      continue;
00839                   }
00840 
00841                   TBranchProxyClassDescriptor *cldesc = 0;
00842 
00843                   if (branchEndname == element->GetName()) {
00844                      // We have a proper node for the base class, recurse
00845 
00846                      if (branch->GetListOfBranches()->GetEntries() == 0) {
00847                         // The branch contains a non-split base class that we are unfolding!
00848 
00849                         // See AnalyzeTree for similar code!
00850                         TBranchProxyClassDescriptor *local_cldesc = 0;
00851 
00852                         TVirtualStreamerInfo *binfo = branch->GetInfo();
00853                         if (strcmp(cl->GetName(),binfo->GetName())!=0) {
00854                            binfo = cl->GetStreamerInfo(); // might be the wrong version
00855                         }
00856                         local_cldesc = new TBranchProxyClassDescriptor(cl->GetName(), binfo,
00857                                                                        branch->GetName(),
00858                                                                        isclones, 0 /* unsplit object */, 
00859                                                                        containerName);
00860 
00861                         TStreamerElement *elem = 0;
00862 
00863                         TIter next(binfo->GetElements());
00864                         while( (elem = (TStreamerElement*)next()) ) {
00865                            AnalyzeElement(branch,elem,level+1,local_cldesc,"");
00866 
00867                         }
00868                         if (NeedToEmulate(cl,0)) {
00869                            proxyTypeName = local_cldesc->GetName();
00870                            local_cldesc = AddClass(local_cldesc);
00871                         }
00872 
00873                      } else {
00874 
00875                         Int_t pos = branchname.Last('.');
00876                         if (pos != -1) {
00877                            branchname.Remove(pos);
00878                         }
00879                         TString local_prefix = topdesc ? topdesc->GetSubBranchPrefix() : parent->GetName();
00880                         cldesc = new TBranchProxyClassDescriptor(cl->GetName(), objInfo,
00881                                                                  branchname,
00882                                                                  local_prefix,
00883                                                                  isclones, branch->GetSplitLevel(),
00884                                                                  containerName);
00885                         lookedAt += AnalyzeBranches( level+1, cldesc, branch, objInfo);
00886                      }
00887                   } else {
00888                      // We do not have a proper node for the base class, we need to loop over
00889                      // the next branches
00890                      Int_t pos = branchname.Last('.');
00891                      if (pos != -1) {
00892                         branchname.Remove(pos);
00893                      }
00894                      TString local_prefix = topdesc ? topdesc->GetSubBranchPrefix() : parent->GetName();
00895                      objInfo = GetBaseClass( element );
00896                      if (objInfo == 0) {
00897                         // There is no data in this base class
00898                         continue;
00899                      }
00900                      cl = objInfo->GetClass();
00901                      cldesc = new TBranchProxyClassDescriptor(cl->GetName(), objInfo,
00902                                                               branchname,
00903                                                               local_prefix,
00904                                                               isclones, branch->GetSplitLevel(),
00905                                                               containerName);
00906                      usedBranch = kFALSE;
00907                      lookedAt += AnalyzeBranches( level, cldesc, branches, objInfo );
00908                   }
00909 
00910                   TBranchProxyClassDescriptor *added = AddClass(cldesc);
00911                   if (added) proxyTypeName = added->GetName();
00912 
00913                } else {
00914                   TBranchProxyClassDescriptor *cldesc = 0;
00915 
00916                   if (branchEndname == element->GetName()) {
00917 
00918                      // We have a proper node for the base class, recurse
00919                      if (branch->GetListOfBranches()->GetEntries() == 0) {
00920                         // The branch contains a non-split object that we are unfolding!
00921 
00922                         // See AnalyzeTree for similar code!
00923                         TBranchProxyClassDescriptor *local_cldesc = 0;
00924 
00925                         TVirtualStreamerInfo *binfo = branch->GetInfo();
00926                         if (strcmp(cl->GetName(),binfo->GetName())!=0) {
00927                            binfo = cl->GetStreamerInfo(); // might be the wrong version
00928                         }
00929                         local_cldesc = new TBranchProxyClassDescriptor(cl->GetName(), binfo,
00930                                                                        branch->GetName(),
00931                                                                        isclones, 0 /* unsplit object */,
00932                                                                        containerName);
00933 
00934                         TStreamerElement *elem = 0;
00935 
00936                         TIter next(binfo->GetElements());
00937                         while( (elem = (TStreamerElement*)next()) ) {
00938                            AnalyzeElement(branch,elem,level+1,local_cldesc,"");
00939                         }
00940 
00941                         if (NeedToEmulate(cl,0)) {
00942                            proxyTypeName = local_cldesc->GetName();
00943                            local_cldesc = AddClass(local_cldesc);
00944                         }
00945 
00946                      } else {
00947 
00948                         if (isclones != TBranchProxyClassDescriptor::kOut) {
00949                            // We have to guess the version number!
00950                            cl = TClass::GetClass(cname);
00951                            objInfo = GetStreamerInfo(branch,branch->GetListOfBranches(),cl);
00952                         }
00953                         cldesc = new TBranchProxyClassDescriptor(cl->GetName(), objInfo,
00954                                                                  branch->GetName(),
00955                                                                  branch->GetName(),
00956                                                                  isclones, branch->GetSplitLevel(),
00957                                                                  containerName);
00958                         lookedAt += AnalyzeBranches( level+1, cldesc, branch, objInfo);
00959                      }
00960                   } else {
00961                      // We do not have a proper node for the base class, we need to loop over
00962                      // the next branches
00963                      TString local_prefix = topdesc ? topdesc->GetSubBranchPrefix() : parent->GetName();
00964                      if (local_prefix.Length()) local_prefix += ".";
00965                      local_prefix += element->GetName();
00966                      objInfo = branch->GetInfo();
00967                      Int_t pos = branchname.Last('.');
00968                      if (pos != -1) {
00969                         branchname.Remove(pos);
00970                      }
00971                      if (isclones != TBranchProxyClassDescriptor::kOut) {
00972                         // We have to guess the version number!
00973                         cl = TClass::GetClass(cname);
00974                         objInfo = GetStreamerInfo(branch, branches, cl);
00975                      }
00976                      cldesc = new TBranchProxyClassDescriptor(cl->GetName(), objInfo,
00977                                                               branchname,
00978                                                               local_prefix,
00979                                                               isclones, branch->GetSplitLevel(),
00980                                                               containerName);
00981                      usedBranch = kFALSE;
00982                      skipped = kTRUE;
00983                      lookedAt += AnalyzeBranches( level, cldesc, branches, objInfo );
00984                   }
00985 
00986                   TBranchProxyClassDescriptor *added = AddClass(cldesc);
00987                   if (added) proxyTypeName = added->GetName();
00988 
00989                }
00990 
00991                AddForward(cl);
00992                AddHeader(cl);
00993                break;
00994             }
00995 
00996             default:
00997                Error("AnalyzeBranch",
00998                      "Unsupported type for %s (%d).",
00999                      branch->GetName(),element->GetType());
01000 
01001          }
01002 
01003          TString dataMemberName = element->GetName();
01004          if (topdesc) {
01005             topdesc->AddDescriptor(  new TBranchProxyDescriptor( dataMemberName.Data(),
01006                proxyTypeName, branchname, true, skipped ), isBase );
01007          } else {
01008             dataMemberName.Prepend(prefix);
01009             AddDescriptor( new TBranchProxyDescriptor( dataMemberName.Data(),
01010                proxyTypeName, branchname, true, skipped ) );
01011          }
01012 
01013          if (usedBranch) {
01014             branches.Next();
01015             ++lookedAt;
01016          }
01017       }
01018       return lookedAt;
01019    }
01020 
01021    UInt_t TTreeProxyGenerator::AnalyzeOldLeaf(TLeaf *leaf, UInt_t /* level */,
01022                                               TBranchProxyClassDescriptor *topdesc)
01023    {
01024       // Analyze the leaf and populate the `TTreeProxyGenerator or
01025       // the topdesc with its findings.
01026 
01027       if (leaf->IsA()==TLeafObject::Class()) {
01028          Error("AnalyzeOldLeaf","TLeafObject not supported yet");
01029          return 0;
01030       }
01031 
01032       TString leafTypeName = leaf->GetTypeName();
01033       Int_t pos = leafTypeName.Last('_');
01034       if (pos!=-1) leafTypeName.Remove(pos);
01035 
01036       Int_t len = leaf->GetLen();
01037       TLeaf *leafcount = leaf->GetLeafCount();
01038 
01039       UInt_t dim = 0;
01040       std::vector<Int_t> maxDim;
01041       //maxDim[0] = maxDim[1] = maxDim[2] = 1;
01042 
01043       TString dimensions;
01044       TString temp = leaf->GetName();
01045       pos = temp.Index("[");
01046       if (pos!=-1) {
01047          if (pos) temp.Remove(0,pos);
01048          dimensions.Append(temp);
01049       }
01050       temp = leaf->GetTitle();
01051       pos = temp.Index("[");
01052       if (pos!=-1) {
01053          if (pos) temp.Remove(0,pos);
01054          dimensions.Append(temp);
01055       }
01056 
01057       Int_t dimlen = dimensions.Length();
01058 
01059       if (dimlen) {
01060          const char *current = dimensions.Data();
01061 
01062          Int_t index;
01063          Int_t scanindex ;
01064          while (current) {
01065             current++;
01066             if (current[0] == ']') {
01067                maxDim.push_back(-1); // maxDim[dim] = -1; // Loop over all elements;
01068             } else {
01069                scanindex = sscanf(current,"%d",&index);
01070                if (scanindex) {
01071                   maxDim.push_back(index); // maxDim[dim] = index;
01072                } else {
01073                   maxDim.push_back(-2); // maxDim[dim] = -2; // Index is calculated via a variable.
01074                }
01075             }
01076             dim ++;
01077             current = (char*)strstr( current, "[" );
01078          }
01079 
01080       }
01081       //char *twodim = (char*)strstr(leaf->GetTitle(),"][");
01082 
01083       if (leafcount) {
01084          len = leafcount->GetMaximum();
01085       }
01086 
01087 
01088       TString type;
01089       switch (dim) {
01090          case 0: {
01091             type = "T";
01092             type += leafTypeName;
01093             type += "Proxy";
01094             break;
01095          }
01096          case 1: {
01097             type = "TArray";
01098             type += leafTypeName;
01099             type += "Proxy";
01100             break;
01101          }
01102          default: {
01103             type = "TArrayProxy<";
01104             for(Int_t ind = dim - 2; ind > 0; --ind) {
01105                type += "TMultiArrayType<";
01106             }
01107             type += "TArrayType<";
01108             type += leaf->GetTypeName();
01109             type += ",";
01110             type += maxDim[dim-1];
01111             type += "> ";
01112             for(Int_t ind = dim - 2; ind > 0; --ind) {
01113                type += ",";
01114                type += maxDim[ind];
01115                type += "> ";
01116             }
01117             type += ">";
01118             break;
01119          }
01120       }
01121 
01122       TString branchName = leaf->GetBranch()->GetName();
01123       TString dataMemberName = leaf->GetName();
01124 
01125       if (topdesc) {
01126          topdesc->AddDescriptor( new TBranchProxyDescriptor( dataMemberName.Data(),
01127                                                              type,
01128                                                              branchName.Data(),
01129                                                              true, false, true ),
01130                                  0 );
01131       } else {
01132          AddDescriptor( new TBranchProxyDescriptor( dataMemberName.Data(),
01133                                                     type,
01134                                                     branchName.Data(),
01135                                                     true, false, true ) );
01136       }
01137 
01138       return 0;
01139 
01140    }
01141 
01142    UInt_t TTreeProxyGenerator::AnalyzeOldBranch(TBranch *branch, UInt_t level,
01143                                                 TBranchProxyClassDescriptor *topdesc)
01144    {
01145       // Analyze the branch and populate the TTreeProxyGenerator or the topdesc with
01146       // its findings.  Sometimes several branch of the mom are also analyzed,
01147       // the number of such branches is returned (this happens in the case of
01148       // embedded objects inside an object inside a clones array split more than
01149       // one level.
01150 
01151       UInt_t extraLookedAt = 0;
01152       TString prefix;
01153 
01154       TString branchName = branch->GetName();
01155 
01156       TObjArray *leaves = branch->GetListOfLeaves();
01157       Int_t nleaves = leaves ? leaves->GetEntriesFast() : 0;
01158 
01159       if (nleaves>1) {
01160 
01161          // Create a holder
01162          TString type = "unknown";
01163          TBranchProxyClassDescriptor *cldesc = AddClass( new TBranchProxyClassDescriptor(branch->GetName()) );
01164          if (cldesc) {
01165             type = cldesc->GetName();
01166 
01167             for(int l=0;l<nleaves;l++) {
01168                TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(l);
01169                extraLookedAt += AnalyzeOldLeaf(leaf,level+1,cldesc);
01170             }
01171          }
01172 
01173          TString dataMemberName = branchName;
01174 
01175          if (topdesc) {
01176             topdesc->AddDescriptor(  new TBranchProxyDescriptor( dataMemberName.Data(),
01177                                                                  type,
01178                                                                  branchName.Data() ),
01179                                     0 );
01180          } else {
01181             // leafname.Prepend(prefix);
01182             AddDescriptor( new TBranchProxyDescriptor( dataMemberName.Data(),
01183                                                        type,
01184                                                        branchName.Data() ) );
01185          }
01186 
01187       } else {
01188 
01189          TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
01190          extraLookedAt += AnalyzeOldLeaf(leaf,level,topdesc);
01191 
01192       }
01193 
01194 
01195       return extraLookedAt;
01196 
01197    }
01198 
01199    void TTreeProxyGenerator::AnalyzeTree(TTree *tree)
01200    {
01201       // Analyze a TTree and its (potential) friends.
01202 
01203       TIter next( tree->GetListOfBranches() );
01204       TBranch *branch;
01205       while ( (branch = (TBranch*)next()) ) {
01206          TVirtualStreamerInfo *info = 0;
01207          const char *branchname = branch->GetName();
01208          const char *classname = branch->GetClassName();
01209          if (classname && strlen(classname)) {
01210             AddForward( classname );
01211             AddHeader( classname );
01212          }
01213 
01214          TBranchProxyClassDescriptor *desc = 0;
01215          TClass *cl = TClass::GetClass(classname);
01216          TString type = "unknown";
01217          if (cl) {
01218             TBranchProxyClassDescriptor::ELocation isclones = TBranchProxyClassDescriptor::kOut;
01219             TString containerName = "";
01220             if (cl==TClonesArray::Class()) {
01221                isclones = TBranchProxyClassDescriptor::kClones;
01222                containerName = "TClonesArray";
01223                if (branch->IsA()==TBranchElement::Class()) {
01224                   const char *cname = ((TBranchElement*)branch)->GetClonesName();
01225                   TClass *ncl = TClass::GetClass(cname);
01226                   if (ncl) {
01227                      cl = ncl;
01228                      info = GetStreamerInfo(branch,branch->GetListOfBranches(),cl);
01229                   } else {
01230                      Error("AnalyzeTree",
01231                            "Introspection of TClonesArray in older file not implemented yet.");
01232                   }
01233                } else {
01234                   TClonesArray **ptr = (TClonesArray**)branch->GetAddress();
01235                   TClonesArray *clones = 0;
01236                   if (ptr==0) {
01237                      clones = new TClonesArray;
01238                      branch->SetAddress(&clones);
01239                      ptr = &clones;
01240                   }
01241                   branch->GetEntry(0);
01242                   TClass *ncl = *ptr ? (*ptr)->GetClass() : 0;
01243                   if (ncl) {
01244                      cl = ncl;
01245                   } else {
01246                      Error("AnalyzeTree",
01247                            "Introspection of TClonesArray for %s failed.",branch->GetName());
01248                   }
01249                }
01250             } else if (cl->GetCollectionProxy()) {
01251                isclones = TBranchProxyClassDescriptor::kSTL;
01252                containerName = cl->GetName();
01253                if (cl->GetCollectionProxy()->GetValueClass()) {
01254                   cl = cl->GetCollectionProxy()->GetValueClass();
01255                } else {
01256                   type = Form("TStlSimpleProxy<%s >", cl->GetName());
01257                   AddHeader(cl);
01258                   if (!cl->IsLoaded()) AddPragma(Form("#pragma link6 C++ class %s;\n", cl->GetName()));
01259                   AddDescriptor( new TBranchProxyDescriptor( branchname, type, branchname ) );
01260                   continue;
01261                }
01262             }
01263             if (cl) {
01264                if (NeedToEmulate(cl,0) || branchname[strlen(branchname)-1] == '.' || branch->GetSplitLevel()) {
01265                   TBranchElement *be = dynamic_cast<TBranchElement*>(branch);
01266                   TVirtualStreamerInfo *beinfo = (be && isclones == TBranchProxyClassDescriptor::kOut) 
01267                      ? be->GetInfo() : cl->GetStreamerInfo(); // the 2nd hand need to be fixed
01268                   desc = new TBranchProxyClassDescriptor(cl->GetName(), beinfo, branchname,
01269                      isclones, branch->GetSplitLevel(),containerName);
01270                } else {
01271                   type = Form("TObjProxy<%s >",cl->GetName());
01272                }
01273             }
01274          }
01275 
01276          if ( branch->GetListOfBranches()->GetEntries() == 0 ) {
01277 
01278             if (cl) {
01279                // We have a non-splitted object!
01280 
01281                if (desc) {
01282                   TVirtualStreamerInfo *cinfo = cl->GetStreamerInfo();
01283                   TStreamerElement *elem = 0;
01284 
01285                   TIter cnext(cinfo->GetElements());
01286                   while( (elem = (TStreamerElement*)cnext()) ) {
01287                      AnalyzeElement(branch,elem,1,desc,"");
01288                   }
01289 
01290                   desc = AddClass(desc);
01291                   if (desc) {
01292                      type = desc->GetName();
01293 
01294                      TString dataMemberName = branchname;
01295 
01296                      AddDescriptor( new TBranchProxyDescriptor( dataMemberName, type, branchname ) );
01297                   }
01298                }
01299             } else {
01300 
01301                // We have a top level raw type.
01302                AnalyzeOldBranch(branch, 0, 0);
01303             }
01304 
01305          } else {
01306 
01307             // We have a splitted object
01308 
01309             TIter subnext( branch->GetListOfBranches() );
01310             if (desc) {
01311                AnalyzeBranches(1,desc,dynamic_cast<TBranchElement*>(branch),info);
01312             }
01313             desc = AddClass(desc);
01314             type = desc->GetName();
01315             TString dataMemberName = branchname;
01316             AddDescriptor( new TBranchProxyDescriptor( dataMemberName, type, branchname ) );
01317 
01318             if ( branchname[strlen(branchname)-1] != '.' ) {
01319                // If there is no dot also include the data member directly
01320 
01321                AnalyzeBranches(1,0,dynamic_cast<TBranchElement*>(branch),info);
01322 
01323                subnext.Reset();
01324             }
01325 
01326          } // if split or non split
01327       }
01328 
01329       // Now let's add the TTreeFriend (if any)
01330       if (tree->GetListOfFriends()) {
01331          TFriendElement *fe;
01332          Int_t count = 0;
01333 
01334          TIter nextfriend(tree->GetListOfFriends());
01335          while ((fe = (TFriendElement*)nextfriend())) {
01336             TTree *t = fe->GetTree();
01337             TFriendProxyDescriptor *desc;
01338             desc = new TFriendProxyDescriptor(t->GetName(), fe->GetName(), count);
01339 
01340             AddFriend( desc );
01341 
01342             fCurrentListOfTopProxies = desc->GetListOfTopProxies();
01343             AnalyzeTree(t);
01344 
01345             count++;
01346          }
01347       }
01348       fCurrentListOfTopProxies = &fListOfTopProxies;
01349    }
01350 
01351    void TTreeProxyGenerator::AnalyzeElement(TBranch *branch, TStreamerElement *element,
01352                                             UInt_t level, TBranchProxyClassDescriptor *topdesc,
01353                                             const char *path)
01354    {
01355       // Analyze the element and populate the TTreeProxyGenerator or the topdesc with
01356       // its findings.
01357 
01358       TString dataMemberName;
01359       TString pxDataMemberName;
01360       TString type;
01361 
01362       // TString prefix;
01363       Bool_t isBase = false;
01364       TString cname;
01365       TString middle;
01366       TBranchProxyClassDescriptor::ELocation isclones = TBranchProxyClassDescriptor::kOut;
01367       TString containerName;
01368       EContainer container = kNone;
01369       if (topdesc) {
01370          if (topdesc->IsClones()) {
01371             container = kClones;
01372             middle = "Cla";
01373             isclones = TBranchProxyClassDescriptor::kClones;
01374             containerName = "TClonesArray";
01375          } else if (topdesc->IsSTL()) {
01376             container = kSTL;
01377             middle = "Stl";
01378             isclones = TBranchProxyClassDescriptor::kSTL;
01379             containerName = topdesc->GetContainerName();
01380          }
01381       }
01382 
01383       if (!element) return;
01384 
01385       if (strcmp(element->GetName(),"This")==0) {
01386          // Skip the artifical streamer element.
01387          return;
01388       }
01389 
01390       if (element->GetType()==-1) {
01391          // This is an ignored TObject base class.
01392          return;
01393       }
01394 
01395 
01396       Bool_t ispointer = false;
01397       switch(element->GetType()) {
01398 
01399          case TVirtualStreamerInfo::kBool:    { type = "T" + middle + "BoolProxy"; break; }
01400          case TVirtualStreamerInfo::kChar:    { type = "T" + middle + "CharProxy"; break; }
01401          case TVirtualStreamerInfo::kShort:   { type = "T" + middle + "ShortProxy"; break; }
01402          case TVirtualStreamerInfo::kInt:     { type = "T" + middle + "IntProxy"; break; }
01403          case TVirtualStreamerInfo::kLong:    { type = "T" + middle + "LongProxy"; break; }
01404          case TVirtualStreamerInfo::kLong64:  { type = "T" + middle + "Long64Proxy"; break; }
01405          case TVirtualStreamerInfo::kFloat:   { type = "T" + middle + "FloatProxy"; break; }
01406          case TVirtualStreamerInfo::kFloat16: { type = "T" + middle + "Float16Proxy"; break; }
01407          case TVirtualStreamerInfo::kDouble:  { type = "T" + middle + "DoubleProxy"; break; }
01408          case TVirtualStreamerInfo::kDouble32:{ type = "T" + middle + "Double32Proxy"; break; }
01409          case TVirtualStreamerInfo::kUChar:   { type = "T" + middle + "UCharProxy"; break; }
01410          case TVirtualStreamerInfo::kUShort:  { type = "T" + middle + "UShortProxy"; break; }
01411          case TVirtualStreamerInfo::kUInt:    { type = "T" + middle + "UIntProxy"; break; }
01412          case TVirtualStreamerInfo::kULong:   { type = "T" + middle + "ULongProxy"; break; }
01413          case TVirtualStreamerInfo::kULong64: { type = "T" + middle + "ULong64Proxy"; break; }
01414          case TVirtualStreamerInfo::kBits:    { type = "T" + middle + "UIntProxy"; break; }
01415 
01416          case TVirtualStreamerInfo::kCharStar: { type = GetArrayType(element,"Char",container); break; }
01417 
01418             // array of basic types  array[8]
01419          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kBool:    { type = GetArrayType(element,"Bool",container ); break; }
01420          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kChar:    { type = GetArrayType(element,"Char",container ); break; }
01421          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kShort:   { type = GetArrayType(element,"Short",container ); break; }
01422          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kInt:     { type = GetArrayType(element,"Int",container ); break; }
01423          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kLong:    { type = GetArrayType(element,"Long",container ); break; }
01424          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kLong64:  { type = GetArrayType(element,"Long64",container ); break; }
01425          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kFloat:   { type = GetArrayType(element,"Float",container ); break; }
01426          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kFloat16: { type = GetArrayType(element,"Float16",container ); break; }
01427          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kDouble:  { type = GetArrayType(element,"Double",container ); break; }
01428          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kDouble32:{ type = GetArrayType(element,"Double32",container ); break; }
01429          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUChar:   { type = GetArrayType(element,"UChar",container ); break; }
01430          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUShort:  { type = GetArrayType(element,"UShort",container ); break; }
01431          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kUInt:    { type = GetArrayType(element,"UInt",container ); break; }
01432          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kULong:   { type = GetArrayType(element,"ULong",container ); break; }
01433          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kULong64: { type = GetArrayType(element,"ULong64",container ); break; }
01434          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kBits:    { type = GetArrayType(element,"UInt",container ); break; }
01435 
01436             // pointer to an array of basic types  array[n]
01437          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kBool:    { type = GetArrayType(element,"Bool",container ); break; }
01438          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kChar:    { type = GetArrayType(element,"Char",container ); break; }
01439          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kShort:   { type = GetArrayType(element,"Short",container ); break; }
01440          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kInt:     { type = GetArrayType(element,"Int",container ); break; }
01441          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kLong:    { type = GetArrayType(element,"Long",container ); break; }
01442          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kLong64:  { type = GetArrayType(element,"Long64",container ); break; }
01443          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kFloat:   { type = GetArrayType(element,"Float",container ); break; }
01444          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kFloat16: { type = GetArrayType(element,"Float16",container ); break; }
01445          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kDouble:  { type = GetArrayType(element,"Double",container ); break; }
01446          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kDouble32:{ type = GetArrayType(element,"Double32",container ); break; }
01447          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUChar:   { type = GetArrayType(element,"UChar",container ); break; }
01448          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUShort:  { type = GetArrayType(element,"UShort",container ); break; }
01449          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kUInt:    { type = GetArrayType(element,"UInt",container ); break; }
01450          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kULong:   { type = GetArrayType(element,"ULong",container ); break; }
01451          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kULong64: { type = GetArrayType(element,"ULong64",container ); break; }
01452          case TVirtualStreamerInfo::kOffsetP + TVirtualStreamerInfo::kBits:    { type = GetArrayType(element,"UInt",container ); break; }
01453 
01454             // array counter //[n]
01455          case TVirtualStreamerInfo::kCounter: { type = "T" + middle + "IntProxy"; break; }
01456 
01457 
01458          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObjectp:
01459          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObjectP:
01460          case TVirtualStreamerInfo::kObjectp:
01461          case TVirtualStreamerInfo::kObjectP:
01462          case TVirtualStreamerInfo::kAnyp:
01463          case TVirtualStreamerInfo::kAnyP:
01464          case TVirtualStreamerInfo::kSTL + TVirtualStreamerInfo::kObjectp:
01465          case TVirtualStreamerInfo::kSTL + TVirtualStreamerInfo::kObjectP:
01466             // set as pointers and fall through to the next switches
01467             ispointer = true;
01468          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kObject:
01469          case TVirtualStreamerInfo::kObject:
01470          case TVirtualStreamerInfo::kTString:
01471          case TVirtualStreamerInfo::kTNamed:
01472          case TVirtualStreamerInfo::kTObject:
01473          case TVirtualStreamerInfo::kAny:
01474          case TVirtualStreamerInfo::kOffsetL + TVirtualStreamerInfo::kAny:
01475          case TVirtualStreamerInfo::kSTL:
01476          case TVirtualStreamerInfo::kBase: {
01477             TClass *cl = element->GetClassPointer();
01478             if (cl) {
01479                type = Form("T%sObjProxy<%s >",
01480                            middle.Data(),cl->GetName());
01481                cname = cl->GetName();
01482                if (cl==TClonesArray::Class()) {
01483                   isclones = TBranchProxyClassDescriptor::kClones;
01484                   containerName = "TClonesArray";
01485 
01486                   Long64_t i = branch->GetTree()->GetReadEntry();
01487                   if (i<0) i = 0;
01488                   branch->GetEntry(i);
01489 
01490                   //char *obj = branch->GetObject();
01491 
01492                   // now need to follow it through to this pointer!
01493 
01494                   TClonesArray *arr;
01495 
01496                   TString fullpath = branch->GetName();
01497                   fullpath += ".";
01498                   if (path && strlen(path)>0) fullpath.Append(path).Append(".");
01499                   fullpath += element->GetName();
01500 
01501                   TTreeFormula *formula = new TTreeFormula("clones",fullpath,branch->GetTree());
01502 
01503                   TFormLeafInfo *leafinfo = formula->GetLeafInfo(0);
01504                   TLeaf *leaf = formula->GetLeaf(0);
01505                   R__ASSERT(leaf && leafinfo);
01506 
01507                   arr = (TClonesArray*)leafinfo->GetLocalValuePointer(leaf,0);
01508 
01509                   /*
01510                     if (ispointer) {
01511                     arr = (TClonesArray*)*(void**)(obj+lOffset);
01512                     } else {
01513                     arr = (TClonesArray*)(obj+lOffset);
01514                     }
01515                   */
01516                   cname = arr->GetClass()->GetName();
01517 
01518                   if (cname.Length()==0) {
01519                      Error("AnalyzeTree",
01520                            "Introspection of TClonesArray in older file not implemented yet.");
01521                   }
01522                   delete formula;
01523                } else if (cl->GetCollectionProxy()) {
01524                   isclones = TBranchProxyClassDescriptor::kSTL;
01525                   containerName = cl->GetName();
01526                   cl = cl->GetCollectionProxy()->GetValueClass();
01527                }
01528             }
01529             else Error("AnalyzeTree","missing class for %s.",branch->GetName());
01530             if (element->IsA()==TStreamerBase::Class()) {
01531                // prefix  = "base";
01532                isBase = true;
01533             }
01534             AddForward(cl);
01535             AddHeader(cl);
01536             break;
01537          }
01538 
01539          default:
01540             Error("AnalyzeTree",
01541                   "Unsupported type for %s %s %d",
01542                   branch->GetName(), element->GetName(), element->GetType());
01543 
01544       }
01545 
01546       dataMemberName = element->GetName();
01547 
01548       if (level<=fMaxUnrolling) {
01549 
01550          // See AnalyzeTree for similar code!
01551          TBranchProxyClassDescriptor *cldesc;
01552 
01553          TClass *cl = TClass::GetClass(cname);
01554          if (cl && cl->CanSplit()) {
01555             cldesc = new TBranchProxyClassDescriptor(cl->GetName(), cl->GetStreamerInfo(),
01556                                                      branch->GetName(),
01557                                                      isclones, 0 /* non-split object */,
01558                                                      containerName);
01559 
01560             TVirtualStreamerInfo *info = cl->GetStreamerInfo();
01561             TStreamerElement *elem = 0;
01562 
01563             TString subpath = path;
01564             if (subpath.Length()>0) subpath += ".";
01565             subpath += dataMemberName;
01566 
01567             TIter next(info->GetElements());
01568             while( (elem = (TStreamerElement*)next()) ) {
01569                AnalyzeElement(branch, elem, level+1, cldesc, subpath.Data());
01570             }
01571 
01572             TBranchProxyClassDescriptor *added = AddClass(cldesc);
01573             if (added) type = added->GetName();
01574          }
01575 
01576       }
01577 
01578       pxDataMemberName = /* prefix + */ dataMemberName;
01579       if (topdesc) {
01580          topdesc->AddDescriptor( new TBranchProxyDescriptor( pxDataMemberName.Data(), type,
01581                                                              dataMemberName.Data(), false),
01582                                  isBase );
01583       } else {
01584          Error("AnalyzeTree","topdesc should not be null in TTreeProxyGenerator::AnalyzeElement.");
01585       }
01586    }
01587 
01588    //----------------------------------------------------------------------------------------------
01589    void TTreeProxyGenerator::ParseOptions()
01590    {
01591       // Parse the options string.
01592 
01593       TString opt = fOptionStr;
01594 
01595       fOptions = 0;
01596       if ( opt.Contains("nohist") ) {
01597          opt.ReplaceAll("nohist","");
01598          fOptions |= kNoHist;
01599       }
01600    }
01601 
01602    //----------------------------------------------------------------------------------------------
01603    static Bool_t R__AddPragmaForClass(TTreeProxyGenerator *gen, TClass *cl)
01604    {
01605       // Add the "pragma C++ class" if needed and return
01606       // true if it has been added _or_ if it is known to
01607       // not be needed.
01608       // (I.e. return kFALSE if a container of this class
01609       // can not have a "pragma C++ class" 
01610       
01611       if (!cl) return kFALSE;
01612       if (cl->GetCollectionProxy()) {
01613          TClass *valcl = cl->GetCollectionProxy()->GetValueClass();
01614          if (!valcl) {
01615             if (!cl->IsLoaded()) gen->AddPragma(Form("#pragma link C++ class %s;\n", cl->GetName()));
01616             return kTRUE;
01617          } else if (R__AddPragmaForClass(gen, valcl)) {
01618             if (!cl->IsLoaded()) gen->AddPragma(Form("#pragma link C++ class %s;\n", cl->GetName()));
01619             return kTRUE;
01620          }
01621       } 
01622       if (cl->IsLoaded()) return kTRUE;
01623       return kFALSE;
01624    }
01625 
01626    //----------------------------------------------------------------------------------------------
01627    static Bool_t R__AddPragmaForClass(TTreeProxyGenerator *gen, const char *classname)
01628    {
01629       // Add the "pragma C++ class" if needed and return
01630       // true if it has been added _or_ if it is known to
01631       // not be needed.
01632       // (I.e. return kFALSE if a container of this class
01633       // can not have a "pragma C++ class" 
01634 
01635       return R__AddPragmaForClass( gen, TClass::GetClass(classname) );
01636 
01637    }
01638 
01639    //----------------------------------------------------------------------------------------------
01640    void TTreeProxyGenerator::WriteProxy()
01641    {
01642       // Check whether the file exist and do something useful if it does
01643       if (fScript.Length()==0) {
01644          Error("WriteProxy","No user script has been specified.");
01645          return;
01646       }
01647 
01648       TString fileLocation = gSystem->DirName(fScript);
01649 
01650       TString incPath = gSystem->GetIncludePath(); // of the form -Idir1  -Idir2 -Idir3
01651       incPath.Append(":").Prepend(" ");
01652       incPath.ReplaceAll(" -I",":");       // of form :dir1 :dir2:dir3
01653       while ( incPath.Index(" :") != -1 ) {
01654          incPath.ReplaceAll(" :",":");
01655       }
01656       incPath.Prepend(fileLocation+":.:");
01657 
01658       const char *filename = gSystem->Which(incPath,fScript);
01659       if (filename==0) {
01660          Error("WriteProxy","Can not find the user's script: %s",fScript.Data());
01661          return;
01662       }
01663       const char *cutfilename = 0;
01664       if (fCutScript.Length()) {
01665          fileLocation = gSystem->DirName(fCutScript);
01666          incPath.Prepend(fileLocation+":.:");
01667          cutfilename = gSystem->Which(incPath,fCutScript);
01668          if (cutfilename==0) {
01669             Error("WriteProxy","Can not find the user's cut script: %s",fCutScript.Data());
01670             delete [] filename;
01671             return;
01672          }
01673       }
01674 
01675       fHeaderFileName = fPrefix;
01676       fHeaderFileName.Append(".h");
01677 
01678       // Check to see if the target file exist.
01679       // If they do we will generate the proxy in temporary file and modify the original
01680       // if and only if it is different.
01681 
01682       Bool_t updating = kFALSE;
01683       if (gSystem->GetPathInfo( fHeaderFileName, 0, (Long_t*)0, 0, 0 ) == 0) {
01684          // file already exist
01685          updating = kTRUE;
01686       }
01687 
01688       TString classname = gSystem->BaseName(fPrefix);
01689 
01690       TString treefile;
01691       Bool_t ischain = fTree->InheritsFrom(TChain::Class());
01692       if (fTree->GetDirectory() && fTree->GetDirectory()->GetFile())
01693          treefile = fTree->GetDirectory()->GetFile()->GetName();
01694       else
01695          treefile = "Memory Directory";
01696 
01697       TString scriptfunc = fScript;
01698       Ssiz_t dot_pos = scriptfunc.Last('.');
01699       if (dot_pos == kNPOS) {
01700          Error("WriteProxy","User's script (%s) has no extension! Nothing will be written.",scriptfunc.Data());
01701          delete [] filename;
01702          delete [] cutfilename;
01703          return;
01704       }
01705       scriptfunc.Replace( dot_pos, fScript.Length()-dot_pos, "");
01706       TString scriptHeader = scriptfunc;
01707       const char * extensions[] = { ".h", ".hh", ".hpp", ".hxx",  ".hPP", ".hXX" };
01708 
01709       int i;
01710       for (i = 0; i < 6; i++ ) {
01711          TString possible = scriptHeader;
01712          possible.Append(extensions[i]);
01713          const char *name = gSystem->Which(incPath,possible);
01714          if (name) {
01715             scriptHeader = possible;
01716             fListOfHeaders.Add(new TNamed("script",Form("#include \"%s\"\n",
01717                                                         scriptHeader.Data())));
01718             delete [] name;
01719             break;
01720          }
01721       }
01722       scriptfunc = gSystem->BaseName(scriptfunc);
01723 
01724 
01725       TString cutscriptfunc = fCutScript;
01726       if (cutfilename) {
01727          dot_pos = cutscriptfunc.Last('.');
01728          cutscriptfunc.Replace( dot_pos, fCutScript.Length()-dot_pos, "");
01729          TString cutscriptHeader = cutscriptfunc;
01730 
01731          for (i = 0; i < 6; i++ ) {
01732             TString possible = cutscriptHeader;
01733             possible.Append(extensions[i]);
01734             const char *name = gSystem->Which(incPath,possible);
01735             if (name) {
01736                cutscriptHeader = possible;
01737                fListOfHeaders.Add(new TNamed("cutscript",Form("#include \"%s\"\n",
01738                                                               cutscriptHeader.Data())));
01739                delete [] name;
01740                break;
01741             }
01742          }
01743          cutscriptfunc = gSystem->BaseName(cutscriptfunc);
01744       }
01745 
01746       FILE *hf;
01747       TString tmpfilename;
01748       if (updating) {
01749          // Coverity[secure_temp]: we don't care about predictable names.
01750          tmpfilename = gSystem->BaseName( tmpnam(0) );
01751          tmpfilename.Append("_proxy.h");
01752          hf = fopen(tmpfilename.Data(),"w");
01753       } else {
01754          hf = fopen(fHeaderFileName.Data(),"w");
01755       }
01756       if (hf == 0) {
01757          Error("WriteProxy","Unable to open the file %s for writing.",fHeaderFileName.Data());
01758          return;
01759       }
01760 
01761       TDatime td;
01762       fprintf(hf,   "/////////////////////////////////////////////////////////////////////////\n");
01763       fprintf(hf,   "//   This class has been automatically generated \n");
01764       fprintf(hf,   "//   (at %s by ROOT version %s)\n",td.AsString(),gROOT->GetVersion());
01765       if (!ischain) {
01766          fprintf(hf,"//   from TTree %s/%s\n",fTree->GetName(),fTree->GetTitle());
01767          fprintf(hf,"//   found on file: %s\n",treefile.Data());
01768       } else {
01769          fprintf(hf,"//   from TChain %s/%s\n",fTree->GetName(),fTree->GetTitle());
01770       }
01771       fprintf(hf,   "/////////////////////////////////////////////////////////////////////////\n");
01772       fprintf(hf,"\n");
01773       fprintf(hf,"\n");
01774 
01775       fprintf(hf,"#ifndef %s_h\n",classname.Data());
01776       fprintf(hf,"#define %s_h\n",classname.Data());
01777       fprintf(hf,"\n");
01778 
01779 
01780       fprintf(hf,"// System Headers needed by the proxy\n");
01781       fprintf(hf,"#if defined(__CINT__) && !defined(__MAKECINT__)\n");
01782       fprintf(hf,"   #define ROOT_Rtypes\n");
01783       fprintf(hf,"   #define ROOT_TError\n");
01784       fprintf(hf,"#endif\n");
01785       fprintf(hf,"#include <TROOT.h>\n");
01786       fprintf(hf,"#include <TChain.h>\n");
01787       fprintf(hf,"#include <TFile.h>\n");
01788       fprintf(hf,"#include <TPad.h>\n");
01789       fprintf(hf,"#include <TH1.h>\n");
01790       fprintf(hf,"#include <TSelector.h>\n");
01791       fprintf(hf,"#include <TBranchProxy.h>\n");
01792       fprintf(hf,"#include <TBranchProxyDirector.h>\n");
01793       fprintf(hf,"#include <TBranchProxyTemplate.h>\n");
01794       fprintf(hf,"#include <TFriendProxy.h>\n");
01795       fprintf(hf,"using namespace ROOT;\n"); // questionable
01796       fprintf(hf,"\n");
01797 
01798       fprintf(hf,"// forward declarations needed by this particular proxy\n");
01799       TIter next( &fListOfForwards );
01800       TObject *current;
01801       while ( (current=next()) ) {
01802          if (strstr(current->GetTitle(),"::")==0) {
01803             // We can not forward declared nested classes (well we might be able to do so for
01804             // the one nested in a namespace but it is not clear yet if we can really reliably
01805             // find this information)
01806             fprintf(hf,"%s",current->GetTitle());
01807          }
01808       }
01809 
01810       fprintf(hf,"\n\n");
01811       fprintf(hf,"// Header needed by this particular proxy\n");
01812       next = &fListOfHeaders;
01813       TObject *header;
01814       while ( (header = next()) ) {
01815          fprintf(hf,"%s",header->GetTitle());
01816       }
01817       fprintf(hf,"\n\n");
01818 
01819       fprintf(hf,"class %s_Interface {\n", scriptfunc.Data());
01820       fprintf(hf,"   // This class defines the list of methods that are directly used by %s,\n",classname.Data());
01821       fprintf(hf,"   // and that can be overloaded in the user's script\n"); 
01822       fprintf(hf,"public:\n");
01823       fprintf(hf,"   void %s_Begin(TTree*) {}\n",scriptfunc.Data());
01824       fprintf(hf,"   void %s_SlaveBegin(TTree*) {}\n",scriptfunc.Data());
01825       fprintf(hf,"   Bool_t %s_Notify() { return kTRUE; }\n",scriptfunc.Data());
01826       fprintf(hf,"   Bool_t %s_Process(Long64_t) { return kTRUE; }\n",scriptfunc.Data());
01827       fprintf(hf,"   void %s_SlaveTerminate() {}\n",scriptfunc.Data());
01828       fprintf(hf,"   void %s_Terminate() {}\n",scriptfunc.Data());
01829       fprintf(hf,"};\n");
01830       fprintf(hf,"\n\n");
01831 
01832       fprintf(hf, "class %s : public TSelector, public %s_Interface {\n", classname.Data(), scriptfunc.Data());
01833       fprintf(hf, "public :\n");
01834       fprintf(hf, "   TTree          *fChain;         //!pointer to the analyzed TTree or TChain\n");
01835       fprintf(hf, "   TH1            *htemp;          //!pointer to the histogram\n");
01836       fprintf(hf, "   TBranchProxyDirector fDirector; //!Manages the proxys\n\n");
01837 
01838       fprintf(hf, "   // Optional User methods\n");
01839       fprintf(hf, "   TClass         *fClass;    // Pointer to this class's description\n");
01840 
01841       if (fListOfClasses.LastIndex()>=0) {
01842          fprintf(hf, "\n   // Wrapper class for each unwounded class\n");
01843          next = &fListOfClasses;
01844          TBranchProxyClassDescriptor *clp;
01845          while ( (clp = (TBranchProxyClassDescriptor*)next()) ) {
01846             clp->OutputDecl(hf, 3, fMaxDatamemberType);
01847          }
01848       }
01849 
01850       if (fListOfFriends.LastIndex()>=0) {
01851          fprintf(hf, "\n   // Wrapper class for each friend TTree\n");
01852          next = &fListOfFriends;
01853          TFriendProxyDescriptor *clp;
01854          while ( (clp = (TFriendProxyDescriptor*)next()) ) {
01855             if (!clp->IsDuplicate()) clp->OutputClassDecl(hf, 3, fMaxDatamemberType);
01856          }
01857       }
01858 
01859       fprintf(hf, "\n   // Proxy for each of the branches, leaves and friends of the tree\n");
01860       next = &fListOfTopProxies;
01861       TBranchProxyDescriptor *data;
01862       while ( (data = (TBranchProxyDescriptor*)next()) ) {
01863          data->OutputDecl(hf, 3, fMaxDatamemberType);
01864       }
01865       if (fListOfFriends.LastIndex()>=0) {
01866          next = &fListOfFriends;
01867          TFriendProxyDescriptor *clp;
01868          while ( (clp = (TFriendProxyDescriptor*)next()) ) {
01869             clp->OutputDecl(hf, 3, fMaxDatamemberType);
01870          }
01871       }
01872       fprintf(hf,"\n\n");
01873 
01874       // Constructor
01875       fprintf(hf,      "   %s(TTree *tree=0) : \n",classname.Data());
01876       fprintf(hf,      "      fChain(0)");
01877       fprintf(hf,   ",\n      htemp(0)");
01878       fprintf(hf,   ",\n      fDirector(tree,-1)");
01879       fprintf(hf,   ",\n      fClass                (TClass::GetClass(\"%s\"))",classname.Data());
01880       next = &fListOfTopProxies;
01881       while ( (data = (TBranchProxyDescriptor*)next()) ) {
01882          fprintf(hf,",\n      %-*s(&fDirector,\"%s\")",
01883                  fMaxDatamemberType, data->GetDataName(), data->GetBranchName());
01884       }
01885       next = &fListOfFriends;
01886       TFriendProxyDescriptor *fpd;
01887       while ( (fpd = (TFriendProxyDescriptor*)next()) ) {
01888           fprintf(hf,",\n      %-*s(&fDirector,tree,%d)",
01889                  fMaxDatamemberType, fpd->GetTitle(), fpd->GetIndex());
01890       }
01891 
01892       fprintf(hf,    "\n      { }\n");
01893 
01894       // Other functions.
01895       fprintf(hf,"   ~%s();\n",classname.Data());
01896       fprintf(hf,"   Int_t   Version() const {return 1;}\n");
01897       fprintf(hf,"   void    Begin(::TTree *tree);\n");
01898       fprintf(hf,"   void    SlaveBegin(::TTree *tree);\n");
01899       fprintf(hf,"   void    Init(::TTree *tree);\n");
01900       fprintf(hf,"   Bool_t  Notify();\n");
01901       fprintf(hf,"   Bool_t  Process(Long64_t entry);\n");
01902       fprintf(hf,"   void    SlaveTerminate();\n");
01903       fprintf(hf,"   void    Terminate();\n");
01904       fprintf(hf,"\n");
01905       fprintf(hf,"   ClassDef(%s,0);\n",classname.Data());
01906       fprintf(hf,"\n\n");
01907 
01908       fprintf(hf,"//inject the user's code\n");
01909       fprintf(hf,"#include \"%s\"\n",fScript.Data());
01910 
01911       if (cutfilename) {
01912          fprintf(hf,"#include \"%s\"\n",fCutScript.Data());
01913       }
01914 
01915       // Close the class.
01916       fprintf(hf,"};\n");
01917       fprintf(hf,"\n");
01918       fprintf(hf,"#endif\n");
01919       fprintf(hf,"\n\n");
01920 
01921       fprintf(hf,"#ifdef __MAKECINT__\n");
01922       if (fListOfClasses.LastIndex()>=0) {
01923          TBranchProxyClassDescriptor *clp;
01924          next = &fListOfClasses;
01925          while ( (clp = (TBranchProxyClassDescriptor*)next()) ) {
01926             fprintf(hf,"#pragma link C++ class %s::%s-;\n",classname.Data(),clp->GetName());
01927             if (clp->GetContainerName().Length()) {
01928                R__AddPragmaForClass(this, clp->GetContainerName());
01929             }
01930          }
01931          next = &fListOfPragmas;
01932          TObjString *prag;
01933          while ( (prag = (TObjString*)next()) ) {
01934             fprintf(hf,"%s",prag->String().Data());
01935          }
01936       }
01937       fprintf(hf,"#pragma link C++ class %s;\n",classname.Data());
01938       fprintf(hf,"#endif\n");
01939       fprintf(hf,"\n\n");
01940 
01941       // Write the implementations.
01942       fprintf(hf,"inline %s::~%s() {\n",classname.Data(),classname.Data());
01943       fprintf(hf,"   // destructor. Clean up helpers.\n");
01944       fprintf(hf,"\n");
01945       fprintf(hf,"}\n");
01946       fprintf(hf,"\n");
01947       fprintf(hf,"inline void %s::Init(TTree *tree)\n",classname.Data());
01948       fprintf(hf,"{\n");
01949       fprintf(hf,"//   Set branch addresses\n");
01950       fprintf(hf,"   if (tree == 0) return;\n");
01951       fprintf(hf,"   fChain = tree;\n");
01952       fprintf(hf,"   fDirector.SetTree(fChain);\n");
01953       fprintf(hf,"   if (htemp == 0) {\n");
01954       fprintf(hf,"      htemp = fDirector.CreateHistogram(GetOption());\n");
01955       if (cutfilename) {
01956          fprintf(hf,"      htemp->SetTitle(\"%s {%s}\");\n",fScript.Data(),fCutScript.Data());
01957       } else {
01958          fprintf(hf,"      htemp->SetTitle(\"%s\");\n",fScript.Data());
01959       }
01960       fprintf(hf,"      fObject = htemp;\n");
01961       fprintf(hf,"   }\n");
01962       fprintf(hf,"}\n");
01963       fprintf(hf,"\n");
01964       fprintf(hf,"Bool_t %s::Notify()\n",classname.Data());
01965       fprintf(hf,"{\n");
01966       fprintf(hf,"   // Called when loading a new file.\n");
01967       fprintf(hf,"   // Get branch pointers.\n");
01968       fprintf(hf,"   fDirector.SetTree(fChain);\n");
01969       fprintf(hf,"   %s_Notify();\n",scriptfunc.Data());
01970       fprintf(hf,"   \n");
01971       fprintf(hf,"   return kTRUE;\n");
01972       fprintf(hf,"}\n");
01973       fprintf(hf,"   \n");
01974 
01975       // generate code for class member function Begin
01976       fprintf(hf,"\n");
01977       fprintf(hf,"inline void %s::Begin(TTree *tree)\n",classname.Data());
01978       fprintf(hf,"{\n");
01979       fprintf(hf,"   // The Begin() function is called at the start of the query.\n");
01980       fprintf(hf,"   // When running with PROOF Begin() is only called on the client.\n");
01981       fprintf(hf,"   // The tree argument is deprecated (on PROOF 0 is passed).\n");
01982       fprintf(hf,"\n");
01983       fprintf(hf,"   TString option = GetOption();\n");
01984       fprintf(hf,"   %s_Begin(tree);\n",scriptfunc.Data());
01985       fprintf(hf,"\n");
01986       fprintf(hf,"}\n");
01987 
01988       // generate code for class member function SlaveBegin
01989       fprintf(hf,"\n");
01990       fprintf(hf,"inline void %s::SlaveBegin(TTree *tree)\n",classname.Data());
01991       fprintf(hf,"{\n");
01992       fprintf(hf,"   // The SlaveBegin() function is called after the Begin() function.\n");
01993       fprintf(hf,"   // When running with PROOF SlaveBegin() is called on each slave server.\n");
01994       fprintf(hf,"   // The tree argument is deprecated (on PROOF 0 is passed).\n");
01995       fprintf(hf,"\n");
01996       fprintf(hf,"   Init(tree);\n");
01997       fprintf(hf,"\n");
01998       fprintf(hf,"   %s_SlaveBegin(tree);\n",scriptfunc.Data());
01999       fprintf(hf,"\n");
02000       fprintf(hf,"}\n");
02001       fprintf(hf,"\n");
02002 
02003       // generate code for class member function Process
02004       fprintf(hf,"inline Bool_t %s::Process(Long64_t entry)\n",classname.Data());
02005       fprintf(hf,"{\n");
02006 
02007       fprintf(hf,"   // The Process() function is called for each entry in the tree (or possibly\n"
02008               "   // keyed object in the case of PROOF) to be processed. The entry argument\n"
02009               "   // specifies which entry in the currently loaded tree is to be processed.\n"
02010               "   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()\n"
02011               "   // to read either all or the required parts of the data. When processing\n"
02012               "   // keyed objects with PROOF, the object is already loaded and is available\n"
02013               "   // via the fObject pointer.\n"
02014               "   //\n"
02015               "   // This function should contain the \"body\" of the analysis. It can contain\n"
02016               "   // simple or elaborate selection criteria, run algorithms on the data\n"
02017               "   // of the event and typically fill histograms.\n\n");
02018       fprintf(hf,"   // WARNING when a selector is used with a TChain, you must use\n");
02019       fprintf(hf,"   //  the pointer to the current TTree to call GetEntry(entry).\n");
02020       fprintf(hf,"   //  The entry is always the local entry number in the current tree.\n");
02021       fprintf(hf,"   //  Assuming that fChain is the pointer to the TChain being processed,\n");
02022       fprintf(hf,"   //  use fChain->GetTree()->GetEntry(entry).\n");
02023       fprintf(hf,"\n");
02024       fprintf(hf,"\n");
02025       fprintf(hf,"   fDirector.SetReadEntry(entry);\n");
02026       if (fOptions & kNoHist) {
02027          if (cutfilename) {
02028             fprintf(hf,"   if (%s()) %s();\n",cutscriptfunc.Data(),scriptfunc.Data());
02029          } else {
02030             fprintf(hf,"   %s();\n",scriptfunc.Data());
02031          }
02032       } else {
02033          if (cutfilename) {
02034             fprintf(hf,"   if (%s()) htemp->Fill(%s());\n",cutscriptfunc.Data(),scriptfunc.Data());
02035          } else {
02036             fprintf(hf,"   htemp->Fill(%s());\n",scriptfunc.Data());
02037          }
02038       }
02039       fprintf(hf,"   %s_Process(entry);\n",scriptfunc.Data());
02040       fprintf(hf,"   return kTRUE;\n");
02041       fprintf(hf,"\n");
02042       fprintf(hf,"}\n\n");
02043 
02044       // generate code for class member function SlaveTerminate
02045       fprintf(hf,"inline void %s::SlaveTerminate()\n",classname.Data());
02046       fprintf(hf,"{\n");
02047       fprintf(hf,"   // The SlaveTerminate() function is called after all entries or objects\n"
02048               "   // have been processed. When running with PROOF SlaveTerminate() is called\n"
02049               "   // on each slave server.");
02050       fprintf(hf,"\n");
02051       fprintf(hf,"   %s_SlaveTerminate();\n",scriptfunc.Data());
02052       fprintf(hf,"}\n\n");
02053 
02054       // generate code for class member function Terminate
02055       fprintf(hf,"inline void %s::Terminate()\n",classname.Data());
02056       fprintf(hf,"{\n");
02057       fprintf(hf,"   // Function called at the end of the event loop.\n");
02058       fprintf(hf,"   htemp = (TH1*)fObject;\n");
02059       fprintf(hf,"   Int_t drawflag = (htemp && htemp->GetEntries()>0);\n");
02060       fprintf(hf,"   \n");
02061       fprintf(hf,"   if (gPad && !drawflag && !fOption.Contains(\"goff\") && !fOption.Contains(\"same\")) {\n");
02062       fprintf(hf,"      gPad->Clear();\n");
02063       fprintf(hf,"   } else {\n");
02064       fprintf(hf,"      if (fOption.Contains(\"goff\")) drawflag = false;\n");
02065       fprintf(hf,"      if (drawflag) htemp->Draw(fOption);\n");
02066       fprintf(hf,"   }\n");
02067       fprintf(hf,"   %s_Terminate();\n",scriptfunc.Data());
02068       fprintf(hf,"}\n");
02069 
02070       fclose(hf);
02071 
02072       if (updating) {
02073          // over-write existing file only if needed.
02074          if (AreDifferent(fHeaderFileName,tmpfilename)) {
02075             gSystem->Unlink(fHeaderFileName);
02076             gSystem->Rename(tmpfilename,fHeaderFileName);
02077          } else gSystem->Unlink(tmpfilename);
02078       }
02079       delete [] filename;
02080       delete [] cutfilename;
02081    }
02082 }

Generated on Tue Jul 5 15:43:13 2011 for ROOT_528-00b_version by  doxygen 1.5.1