// // Program to convert rootples to ntuples and viceversa for the PAMELA ground data. // Written by Emiliano Mocchiutti. // // Version 3.00 (2005/11/29) // // Changelog: // // 2.00 - 3.00 (2005/11/29): compiled. // // 1.05 - 2.00 (2005/10/07): added TOF and TRIGGER level1 conversion (from PAW to ROOT). // // 1.04 - 1.05 (2005/10/04): tracker version 2.00 conversion. // // 1.03 - 1.04 (2005/09/07): small bug unloading libraries fixed. // // 1.02 - 1.03 (2005/08/03): changes for working on 64 bit machines. // // 1.01 - 1.02 (2005/07/21): don't load yodaUtility.c anymore and use clone routines in CaloFunctions.h // // 1.00 - 1.01 (2005/07/14): small change in the call to getLEVname (changed to be compiled). // // 0.1 - 1.00 (2005/07/05): working, it converts AC level1 rootples to ntuples and TRK level1 and level2 ntuples to rootples. // // v. 0.1 (2005/06/15): created. // // #include #include // #include #include #include #include // //#include #include #include #include // #if !defined (__CINT__) #include #include #include #include #include #include extern const char *pathtocalibration(); extern void coptrklev2(const char [], struct Tracklev2 &, int &); extern void cretrklev2(int &, Tracklev2 &); extern void ccltrklev2(struct Tracklev2 &); extern void coptrklev1(const char [], struct Tracklev1 &, int &); extern void cretrklev1(int &, Tracklev1 &); extern void ccltrklev1(struct Tracklev1 &); extern void coptoflev1(const char [], struct Toflev1 &, int &); extern void cretoflev1(int &, Toflev1 &); extern void ccltoflev1(struct Toflev1 &); extern void copaclev1(const char []); extern void cfiaclev1(struct AClev1 &); extern void cclaclev1(); #endif // // #include // #include #include #include int GroundDataConvert(TString filename, TString detector, TString level, TString outDir = "", Int_t FORCE = 0){ gROOT->GetListOfCanvases()->Delete(); gDirectory->GetList()->Delete(); #if defined (__CINT__) emicheckLib(); #endif const char* startingdir = gSystem->WorkingDirectory(); TString path; stringcopy(path,startingdir); const char *decname = detector; TString level2; stringcopy(level2,level); const char *lev = level2; Int_t lev1 = 0; Int_t lev2 = 0; if ( !strcmp(lev,"1") ){ lev1 = 1; }; if ( !strcmp(lev,"2") ){ lev2 = 1; }; Int_t TRK = 0; stringstream bdir; bdir.str(""); if ( !strcmp(decname,"tracker") ) { #if defined (__CINT__) const char *sdir=gSystem->Getenv("PAM_LIB"); if ( !strcmp(lev,"1") ){ bdir << sdir << "/liboptrklev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptrklev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); }; if ( !strcmp(lev,"2") ){ bdir << sdir << "/liboptrklev2.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptrklev2_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev2.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev2_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev2.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev2_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); }; #endif TRK = 1; }; Int_t AC = 0; if ( !strcmp(decname,"anticounter") ) { if ( strcmp(lev,"1") ) printf("WARNING: this program can convert only LEVEL1 AC rootples!\nWARNING: assuming you are meaning level=1\n"); lev="1"; lev1 = 1; lev2 = 0; level = "1"; #if defined (__CINT__) const char *sdir=gSystem->Getenv("PAM_LIB"); bdir << sdir << "/libopaclev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libopaclev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libfiaclev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libfiaclev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libclaclev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libclaclev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); #endif AC = 1; }; Int_t TOF = 0; if ( !strcmp(decname,"tof") ) { if ( strcmp(lev,"1") ) printf("WARNING: this program can convert only LEVEL1 TOF ntuples!\nWARNING: assuming you are meaning level=1\n"); lev="1"; lev1 = 1; lev2 = 0; level = "1"; #if defined (__CINT__) const char *sdir=gSystem->Getenv("PAM_LIB"); bdir << sdir << "/liboptoflev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptoflev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretoflev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretoflev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltoflev1.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltoflev1_C.so"; gSystem->Load(bdir.str().c_str()); bdir.str(""); #endif TOF = 1; }; // printf("\n Check the existence of input data: \n\n"); // Int_t isfileinput = 0; stringstream ifile; if ( TRK ){ const string fil = (const char*)filename; // Int_t posiz = fil.find(".rz"); if ( posiz == -1 ) posiz = fil.find(".RZ"); if ( posiz == -1 ){ const string fil2 = (const char*)filename; Int_t nposiz = fil2.find("dw_"); if ( nposiz == -1 ) nposiz = fil2.find("DW_"); if ( nposiz == -1 ) return(5); nposiz = nposiz+3; Int_t nposiz2 = nposiz+10; TString ufile2; stringcopy(ufile2,filename,nposiz,nposiz2); const char *trname = ufile2; // ifile << filename.Data() << "/Physics/Level"; ifile << lev << "/DW_"; ifile << trname << "_level"; ifile << lev << ".rz"; } else { if ( outDir == "" ){ printf("ERROR: if you give as input the name of the file to convert you must give also the output directory where to store data\n"); printf("Exiting...\n"); return(6); }; TString ufile2; Int_t nposiz = posiz -1; Int_t nposiz2 = posiz; stringcopy(ufile2,filename,nposiz,nposiz2); const char *levname = ufile2; if ( strcmp(lev,levname) ){ printf("ERROR: you asked to process LEVEL%s data but the input file contains LEVEL%s data!\n\n Exiting... \n",lev,levname); return(6); }; ifile << filename; isfileinput = 1; }; ifstream mypawfile; mypawfile.open(ifile.str().c_str()); if ( !mypawfile ){ printf("ERROR: file %s does not exist \n\n",ifile.str().c_str()); printf("Exiting...\n"); return(2); }; printf(" OK, TRACKER input file exists \n\n"); }; // if ( TOF ){ const string fil = (const char*)filename; // Int_t posiz = fil.find(".rz"); if ( posiz == -1 ) posiz = fil.find(".RZ"); if ( posiz == -1 ){ const string fil2 = (const char*)filename; Int_t nposiz = fil2.find("dw_"); if ( nposiz == -1 ) nposiz = fil2.find("DW_"); if ( nposiz == -1 ) return(5); nposiz = nposiz+3; Int_t nposiz2 = nposiz+10; TString ufile2; stringcopy(ufile2,filename,nposiz,nposiz2); const char *trname = ufile2; // ifile << filename.Data() << "/Physics/Level"; ifile << lev << "/DW_"; ifile << trname << "_toftrack.rz"; } else { if ( outDir == "" ){ printf("ERROR: if you give as input the name of the file to convert you must give also the output directory where to store data\n"); printf("Exiting...\n"); return(6); }; // ifile << filename; isfileinput = 1; }; ifstream mypawfile; mypawfile.open(ifile.str().c_str()); if ( !mypawfile ){ printf("ERROR: file %s does not exist \n\n",ifile.str().c_str()); printf("Exiting...\n"); return(2); }; printf(" OK, TOF input file exists \n\n"); }; // TFile *headerFile = 0; TTree *tr = 0; TFile *acFile = 0; if ( AC ){ headerFile = emigetFile(filename, "Physics", "Header"); if ( !headerFile ){ headerFile->Close(); printf("No header file, exiting...\n"); return(1); }; tr = (TTree*)headerFile->Get("Pscu"); // acFile = emigetFile(filename, "Physics.Level1","Anticounter"); if ( !acFile ){ printf("ERROR: No AC file! \n Exiting... \n\n"); return(2); }; printf(" OK, AC file exists \n\n"); tr->AddFriend("AcLevel1",acFile); }; // if ( !FORCE ) printf(" Not in FORCE mode, check the existence of output data: \n"); Int_t nofile = 0; // TString filety = ""; TString detc = ""; // detector // if ( TRK ){ filety = "root"; detc = "Tracker"; // detector }; if ( TOF ){ filety = "root"; detc = "TofTrigger"; // detector }; if ( AC ){ filety = "rz"; detc = "Anticounter"; // detector } char *file; if ( !isfileinput ){ file = getLEVname(filename,detc,lev,filety); } else { file = getfileLEVname(filename,detc,lev,filety); }; // stringstream file2; file2.str(""); stringstream file3; file3.str(""); const char *file4 = 0; if ( outDir == "" ){ file4 = filename; file3 << file4 << "/Physics/Level"; file3 << lev; } else { file4 = outDir; file3 << file4; }; file2 << file3.str().c_str() << "/"; file2 << file; // printf("\n Filename will be: \n %s \n\n",file2.str().c_str()); ifstream mypawfile; mypawfile.open(file2.str().c_str()); if ( mypawfile ){ nofile = 1; } else { if ( !FORCE ) printf("Error in opening file: file %s does not exist \n",file2.str().c_str()); }; // if ( !FORCE ){ if ( nofile ){ printf(" ERROR: file already exists! Use FORCE = 1 to override \n\n"); gSystem->ChangeDirectory(path); return(3); } else { printf("\n OK, I will create it!\n\n"); }; }; // // // if ( TRK ){ if ( lev1 ){ // // create output file: // struct Tracklev1 trklev1; char *type; type = "NEW"; if ( FORCE ) type = "RECREATE"; TFile *hfile = 0; TTree *tree = 0; hfile = new TFile(file2.str().c_str(),type,"Tracker LEVEL1 data"); tree = new TTree("TrkLevel1","PAMELA level1 tracker data"); tree->Branch("good1",&trklev1.good1,"good1/O"); tree->Branch("nev1",&trklev1.nev1,"nev1/I"); tree->Branch("pkt_type1",&trklev1.pkt_type1,"pkt_type1/I"); tree->Branch("pkt_num1",&trklev1.pkt_num1,"pkt_num1/I"); tree->Branch("obt1",&trklev1.obt1,"obt1/I"); tree->Branch("which_calib1",&trklev1.which_calib1,"which_calib1/I"); tree->Branch("nclstr1",&trklev1.nclstr1,"nclstr1/I"); tree->Branch("view",trklev1.view,"view[nclstr1]/I"); tree->Branch("ladder",trklev1.ladder,"ladder[nclstr1]/I"); tree->Branch("maxs",trklev1.maxs,"maxs[nclstr1]/I"); tree->Branch("mult",trklev1.mult,"mult[nclstr1]/I"); tree->Branch("dedx",trklev1.dedx,"dedx[nclstr1]/F"); tree->Branch("indstart",trklev1.indstart,"indstart[nclstr1]/I"); tree->Branch("indmax",trklev1.indmax,"indmax[nclstr1]/I"); tree->Branch("totcllength",&trklev1.totcllength,"totcllength/I"); tree->Branch("clsignal",trklev1.clsignal,"clsignal[totcllength]/F"); tree->Branch("cnev",trklev1.cnev,"cnev[24][12]/F"); // Int_t trnev = 0; coptrklev1(ifile.str().c_str(),trklev1,trnev); printf(" The tracker ntuple contains %i events\n Processing data: \n",trnev ); for ( Int_t i = 1; i < trnev+1; i++){ if ( i%1000 == 0 ) printf("%iK \n",i/1000); // cretrklev1(i,trklev1); // tree->Fill(); }; printf(" Finished! Closing ntuple...\n"); ccltrklev1(trklev1); // hfile->Write(); hfile->Close(); // gSystem->ChangeDirectory(path); #if defined (__CINT__) char *sdir=gSystem->Getenv("PAM_LIB"); bdir.str(""); bdir << sdir << "/libcltrklev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptrklev1_C.so"; gSystem->Unload(bdir.str().c_str()); stringstream cdir; cdir.str(""); char *sdir=gSystem->Getenv("PAM_LIB"); cdir << sdir << "/liboptrklev1.so"; gSystem->Unload(cdir.str().c_str()); #endif return(0); }; if ( lev2 ){ struct Tracklev2 trklev2; char *type; type = "NEW"; if ( FORCE ) type = "RECREATE"; TFile *hfile = 0; TTree *tree = 0; hfile = new TFile(file2.str().c_str(),type,"Tracker LEVEL2 data"); tree = new TTree("TrkLevel2","PAMELA level2 tracker data"); tree->Branch("good2",&trklev2.good2,"good2/O"); tree->Branch("nev2",&trklev2.nev2,"nev2/I"); tree->Branch("pkt_type",&trklev2.pkt_type,"pkt_type/I"); tree->Branch("pkt_num",&trklev2.pkt_num,"pkt_num/I"); tree->Branch("obt",&trklev2.obt,"obt/I"); tree->Branch("which_calib",&trklev2.which_calib,"which_calib/I"); tree->Branch("ntrk",&trklev2.ntrk,"ntrk/I"); tree->Branch("image",&trklev2.image,"image[ntrk]/I"); tree->Branch("xm",trklev2.xm,"xm[ntrk][6]/F"); tree->Branch("ym",trklev2.ym,"ym[ntrk][6]/F"); tree->Branch("zm",trklev2.zm,"zm[ntrk][6]/F"); tree->Branch("resx",trklev2.resx,"resx[ntrk][6]/F"); tree->Branch("resy",trklev2.resy,"resy[ntrk][6]/F"); tree->Branch("al",trklev2.al,"al[ntrk][5]/F"); tree->Branch("coval",trklev2.coval,"coval[ntrk][5][5]/F"); tree->Branch("chi2",trklev2.chi2,"chi2[ntrk]/F"); tree->Branch("xgood",trklev2.xgood,"xgood[ntrk][6]/I"); tree->Branch("ygood",trklev2.ygood,"ygood[ntrk][6]/I"); tree->Branch("xv",trklev2.xv,"xv[ntrk][6]/F"); tree->Branch("yv",trklev2.yv,"yv[ntrk][6]/F"); tree->Branch("zv",trklev2.zv,"zv[ntrk][6]/F"); tree->Branch("axv",trklev2.axv,"axv[ntrk]/F"); tree->Branch("ayv",trklev2.ayv,"ayv[ntrk][6]/F"); tree->Branch("dedxp",trklev2.dedxp,"dedxp[ntrk][6]/F"); tree->Branch("nclsx",trklev2.nclsx,"nclsx[6]/I"); tree->Branch("nclsy",trklev2.nclsy,"nclsy[6]/I"); // Int_t trnev = 0; coptrklev2(ifile.str().c_str(),trklev2,trnev); printf(" The tracker ntuple contains %i events\n Processing data: \n",trnev ); for ( Int_t i = 1; i < trnev+1; i++){ if ( i%1000 == 0 ) printf("%iK \n",i/1000); // cretrklev2(i,trklev2); // tree->Fill(); }; printf(" Finished! Closing ntuple...\n"); ccltrklev2(trklev2); // hfile->Write(); hfile->Close(); // gSystem->ChangeDirectory(path); #if defined (__CINT__) char *sdir=gSystem->Getenv("PAM_LIB"); bdir.str(""); bdir << sdir << "/libcltrklev2_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltrklev2.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev2_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretrklev2.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptrklev2_C.so"; gSystem->Unload(bdir.str().c_str()); stringstream bdir; bdir.str(""); char *sdir=gSystem->Getenv("PAM_LIB"); bdir << sdir << "/liboptrklev2.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); #endif return(0); }; }; if ( TOF ){ if ( lev1 ){ // // create output file: // struct Toflev1 toflev1; char *type; type = "NEW"; if ( FORCE ) type = "RECREATE"; TFile *hfile = 0; TTree *tree = 0; hfile = new TFile(file2.str().c_str(),type,"Tof and Trigger LEVEL1 data"); tree = new TTree("TofLevel1","PAMELA level1 tof and trigger data"); tree->Branch("good2",&toflev1.good2,"good2/O"); tree->Branch("nev2",&toflev1.nev2,"nev2/I"); tree->Branch("pkt_type",&toflev1.pkt_type,"pkt_type/I"); tree->Branch("pkt_num",&toflev1.pkt_num,"pkt_num/I"); tree->Branch("obt",&toflev1.obt,"obt/I"); tree->Branch("which_calib",&toflev1.which_calib,"which_calib/I"); tree->Branch("tdcid",toflev1.tdcid,"tdcid[12]/I"); tree->Branch("evcount",toflev1.evcount,"evcount[12]/I"); tree->Branch("tdcmask",toflev1.tdcmask,"tdcmask[12]/I"); tree->Branch("adc",toflev1.adc,"adc[12][4]/I"); tree->Branch("tdc",toflev1.tdc,"tdc[12][4]/I"); tree->Branch("temp1",toflev1.temp1,"temp1[12]/I"); tree->Branch("temp2",toflev1.temp2,"temp2[12]/I"); tree->Branch("beta",toflev1.beta,"beta[5]/F"); tree->Branch("xtof",toflev1.xtof,"xtof[3]/F"); tree->Branch("ytof",toflev1.ytof,"ytof[3]/F"); tree->Branch("adc_c",toflev1.adc_c,"adc_c[12][4]/F"); tree->Branch("iflag",toflev1.iflag,"iflag[6]/I"); tree->Branch("jflag",toflev1.jflag,"jflag[6]/I"); tree->Branch("xout",toflev1.xout,"xout[3]/F"); tree->Branch("yout",toflev1.yout,"yout[3]/F"); tree->Branch("trig_evcount",&toflev1.trig_evcount,"trig_evcount/I"); tree->Branch("pmtpl",toflev1.pmtpl,"pmtpl[3]/I"); tree->Branch("trigrate",toflev1.trigrate,"trigrate[6]/I"); tree->Branch("dltime",toflev1.dltime,"dltime[2]/I"); tree->Branch("s4calcount",toflev1.s4calcount,"s4calcount[2]/I"); tree->Branch("pmtcount1",toflev1.pmtcount1,"pmtcount1[24]/I"); tree->Branch("pmtcount2",toflev1.pmtcount2,"pmtcount2[24]/I"); tree->Branch("patternbusy",toflev1.patternbusy,"patternbusy[3]/I"); tree->Branch("patterntrig",toflev1.patterntrig,"patterntrig[6]/I"); tree->Branch("trigconf",&toflev1.trigconf,"trigconf/I"); tree->Branch("ntrk",&toflev1.ntrk,"ntrk/I"); tree->Branch("image",&toflev1.image,"image[ntrk]/I"); tree->Branch("xm",toflev1.xm,"xm[ntrk][6]/F"); tree->Branch("ym",toflev1.ym,"ym[ntrk][6]/F"); tree->Branch("zm",toflev1.zm,"zm[ntrk][6]/F"); tree->Branch("resx",toflev1.resx,"resx[ntrk][6]/F"); tree->Branch("resy",toflev1.resy,"resy[ntrk][6]/F"); tree->Branch("al",toflev1.al,"al[ntrk][5]/F"); tree->Branch("coval",toflev1.coval,"coval[ntrk][5][5]/F"); tree->Branch("chi2",toflev1.chi2,"chi2[ntrk]/F"); tree->Branch("xgood",toflev1.xgood,"xgood[ntrk][6]/I"); tree->Branch("ygood",toflev1.ygood,"ygood[ntrk][6]/I"); tree->Branch("xv",toflev1.xv,"xv[ntrk][6]/F"); tree->Branch("yv",toflev1.yv,"yv[ntrk][6]/F"); tree->Branch("zv",toflev1.zv,"zv[ntrk][6]/F"); tree->Branch("axv",toflev1.axv,"axv[ntrk]/F"); tree->Branch("ayv",toflev1.ayv,"ayv[ntrk][6]/F"); tree->Branch("dedxp",toflev1.dedxp,"dedxp[ntrk][6]/F"); tree->Branch("nclsx",toflev1.nclsx,"nclsx[6]/I"); tree->Branch("nclsy",toflev1.nclsy,"nclsy[6]/I"); // Int_t trnev = 0; coptoflev1(ifile.str().c_str(),toflev1,trnev); printf(" The tof ntuple contains %i events\n Processing data: \n",trnev ); for ( Int_t i = 1; i < trnev+1; i++){ if ( i%1000 == 0 ) printf("%iK \n",i/1000); // cretoflev1(i,toflev1); // tree->Fill(); }; printf(" Finished! Closing ntuple...\n"); ccltoflev1(toflev1); // hfile->Write(); hfile->Close(); // gSystem->ChangeDirectory(path); #if defined (__CINT__) char *sdir=gSystem->Getenv("PAM_LIB"); bdir.str(""); bdir << sdir << "/libcltoflev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libcltoflev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretoflev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libretoflev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/liboptoflev1_C.so"; gSystem->Unload(bdir.str().c_str()); stringstream cdir; cdir.str(""); char *sdir=gSystem->Getenv("PAM_LIB"); cdir << sdir << "/liboptoflev1.so"; gSystem->Unload(cdir.str().c_str()); #endif return(0); }; }; if ( AC ){ Long64_t nevents = tr->GetEntries(); if ( nevents < 1 ) { printf("The file is empty!\n"); return(1); }; AnticounterLevel1 *ac = new AnticounterLevel1(); tr->SetBranchAddress("Event",&ac); // // Open the rz file with a fortran call // struct AClev1 aclev1; copaclev1(file2.str().c_str()); // // run over all the events: // for (Int_t i = 0; i < nevents; i++){ if ( i%1000 == 0 ) printf("%iK \n",i/1000); // // open the rootple and copy variables values in "evento_" // tr->GetEntry(i); // // here you must run over all your variables and copy them in evento_.... // aclev1.obt = ac->obt; aclev1.evfile = ac->evfile; aclev1.headc = ac->headc; for (Int_t j = 0; j<2 ; j++){ aclev1.status[j] = (int)ac->status[j]; aclev1.hitmap[j] = (int)ac->hitmap[j]; aclev1.hitstatus[j] = (int)ac->hitstatus[j]; aclev1.trigger[j] = (int)ac->trigger[j]; }; // // fill the ntuple // cfiaclev1(aclev1); }; // // close the rz file // printf("\nClose rz file! \n\n"); cclaclev1(); printf("File saved in \n\n %s \n\n",file2.str().c_str()); gSystem->ChangeDirectory(path); #if defined (__CINT__) char *sdir=gSystem->Getenv("PAM_LIB"); bdir.str(""); bdir << sdir << "/libclaclev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libclaclev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libfiaclev1_C.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libfiaclev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); bdir << sdir << "/libopaclev1_C.so"; gSystem->Unload(bdir.str().c_str()); char *sdir=gSystem->Getenv("PAM_LIB"); stringstream bdir; bdir.str(""); bdir << sdir << "/libopaclev1.so"; gSystem->Unload(bdir.str().c_str()); bdir.str(""); #endif return(0); }; return(0); }