Feature | Status | NOTES |
Class to communicate with R | DONE | class TRInterface |
Class to pass R objects to ROOT objects | DONE | class TRObjectProxy |
Cint support | DONE | You can call it from ROOT's interpreter |
Linux | DONE | Run under serveral linux flavors |
MacOSX | DONE | Tested compilation in macosx(Snow Leopard) |
Windows | PARTIALLY | Tested in cygwin |
ACLiC | PARTIALLY | Using custom rootlogon.C |
Passing of functions from ROOT to R | PARTIALLY | Require ACLiC and TF1,TF2,TF3 and TFormula is not supported |
Overload operators for facilities | PARTIALLY | You can use operator [] for assignation |
cmake | PARTIALLY | Needs to be improved to use R commands to get flags |
testsuite | PROGRESS | All features are being tested with ROOT scripts and compiled code |
#apt-get install r-base r-base-dev
#port install R
http://cygwin.com/setup-x86.exe
gcc gcc-g++ gcc-fortran make xterm git libX11-devel libXft-devel libXext-devel libXpm-deve libXt-deve libGL-devel libGLU-devel X-start-menu-icons xorg-scripts xorg-xserver xorg-server-devel gsl-devel cygutils-x11 libiconv R
R -e "install.packages(\"Rcpp_0.10.4-patched.tar.gz\")"
R -e "install.packages(\"r/src/Rcpp_0.10.6.tar.gz\")" R -e "install.packages(\"r/src/RInside_0.2.10.tar.gz\")"
git clone https://github.com/omazapa/root.git git checkout v5-34-00-patches-r ./configure --enable-r make
./configure win32gcc --enable-r
source $ROOTSYS/bin/thisroot.sh
#include<TRInterface.h> gR->Parse("print(seq(1,5,0.5))"); Int_t power=gR->ParseEval("2^3").ToScalar(); cout<<power<<endl; TVector v=gR->ParseEval("seq(1,5,0.5)").ToVector(); v.Print()
#include<TRInterface.h> ROOT::R::TRInterface &r=gR->Instance(); //non-pointer instance to use operators with a nice syntax. Int_t n=2; TVectorD v(n); v[0]=0; v[1]=1; TMatrixD m(n,n); m[0]=v; m[1]=v; r.Assign(n,"n");//calling the method r["v"]=v; //using operators (*gR)["m"]=m; //using operators with global pointer object.
#include<TRInterface.h> ROOT::R::TRObjectProxy p=gR->ParseEval("2^3"); Int_t power=p.ToScalar();//Gettting scalar type from TVectorD v=gR->ParseEval("seq(1,5,0.5)").ToVector();//getting directly the vector without an TRObjectProxy's object cout<<power<<endl; v.Print();
{ #include <compilerdata.h> gSystem->AddIncludePath(RINCLUDEPATH); gSystem->AddLinkedLibs(RLINKEDLIBS); }
//numerical integration using R //passing the function from ROOT #include<TMath.h> #include<TRInterface.h> #include<Math/Integrator.h> #include<TF1.h> //To integrate using R the function must be vectorized. //The idea is to receive a vector as an argument to evaluate //every element, saving the result in other vector //and return the resultant vector. TVectorD BreitWignerVectorized(TVectorD xx) { TVectorD result(xx.GetNoElements()); for(Int_t i=0;i<xx.GetNoElements();i++) { result[i]=TMath::BreitWigner(xx[i]); } return result; } double BreitWignerWrap( double x){ return TMath::BreitWigner(x); } void Integration() { #if defined(__CINT__) && !defined(__MAKECINT__) cout << "WARNING: This tutorial can run only using ACliC, you must run it by doing: " << endl; cout << "cd $ROOTSYS/tutorials/r/" << endl; cout << "\t .x Integration.C+" << endl; return; #endif ROOT::R::TRInterface &r=gR->Instance(); r["BreitWigner"]=ROOT::R::TRFunction(BreitWignerVectorized); Double_t value=r.ParseEval("integrate(BreitWigner, lower = -2, upper = 2)$value"); std::cout.precision(18); std::cout<<"Integral of BreitWigner Function in the interval [-2, 2] R = "<<value<<std::endl; ROOT::Math::WrappedFunction<> wf(BreitWignerWrap); ROOT::Math::Integrator i(wf); value=i.Integral(-2,2); std::cout<<"Integral of BreitWigner Function in the interval [-2, 2] MathMore = "<<value<<std::endl; TF1 f1("BreitWigner","BreitWignerWrap(x)"); value=f1.Integral(-2,2); std::cout<<"Integral of BreitWigner Function in the interval [-2, 2] TF1 = "<<value<<std::endl; //infinte limits value=r.ParseEval("integrate(BreitWigner, lower = -Inf, upper = Inf)$value"); std::cout<<"Integral of BreitWigner Function in the interval [-Inf, Inf] R = "<<value<<std::endl; }
root -l Integration.C+
[omazapa] [tuxhome] [~]$ root -l root [0] #include<TRInterface.h> root [1] gR->Interactive() [r]:a=seq seq seq_along seq.Date seq.default seq.int seq_len seq.POSIXt sequence [r]:a=seq(1,5,0.5) [r]:.q root [2] TVectorD v=gR->ParseEval("a").ToVector(); root [3] v.Print() Vector (9) is as follows | 1 | ------------------ 0 |1 1 |1.5 2 |2 3 |2.5 4 |3 5 |3.5 6 |4 7 |4.5 8 |5 root [4]
#include<TRInterface.h> TCanvas *Fitting(){ TCanvas *c1 = new TCanvas("c1","Curve Fit",700,500); c1->SetGrid(); // draw a frame for multiples graphs TMultiGraph *mg = new TMultiGraph(); // create the first graph (points with gaussian noise) const Int_t n = 24; Double_t x1[n] ; Double_t y1[n] ; //Generate points along a X^3 with noise TRandom rg; rg.SetSeed(520); for (Int_t i = 0; i < n; i++) { x1[i] = rg.Uniform(0, 1); y1[i] = TMath::Power(x1[i], 3) + rg.Gaus() * 0.06; } TGraph *gr1 = new TGraph(n,x1,y1); gr1->SetMarkerColor(kBlue); gr1->SetMarkerStyle(8); gr1->SetMarkerSize(1); mg->Add(gr1); // create second graph TF1 *f_known=new TF1("f_known","pow(x,3)",0,1); TGraph *gr2 = new TGraph(f_known); gr2->SetMarkerColor(kRed); gr2->SetMarkerStyle(8); gr2->SetMarkerSize(1); mg->Add(gr2); //passing x and y values to R for fitting gR->Assign(TVectorD(n, x1), "x"); gR->Assign(TVectorD(n, y1), "y"); //creating a R data frame gR->Parse("ds<-data.frame(x=x,y=y)"); //fitting x and y to X^power using Nonlinear Least Squares gR->Parse("m <- nls(y ~ I(x^power),data = ds, start = list(power = 1),trace = T)"); //getting the fitted value (power) ROOT::R::TRObjectProxy robj=gR->ParseEval("summary(m)$coefficients[1]"); Double_t power=robj.ToScalar(); TF1 *f_fitted=new TF1("f_fitted","pow(x,[0])",0,1); f_fitted->SetParameter(0,power); //plotting the fitted function TGraph *gr3 = new TGraph(f_fitted); gr3->SetMarkerColor(kGreen); gr3->SetMarkerStyle(8); gr3->SetMarkerSize(1); mg->Add(gr3); mg->Draw("ap"); //displaying basic results TPaveText *pt = new TPaveText(0.1,0.6,0.5,0.9,"brNDC"); pt->SetFillColor(18); pt->SetTextAlign(12); pt->AddText("Fitting x^power "); pt->AddText(" \"Blue\" Points with gaussian noise to be fitted"); pt->AddText(" \"Red\" Known function x^3"); TString fmsg; fmsg.Form(" \"Green\" Fitted function with power=%.4lf",power); pt->AddText(fmsg); pt->Draw(); c1->Update(); return c1; }
root -l Fitting.C
#include<TRInterface.h> //in the next function the pointer *double should be change by TVectorD, because the pointer has no //sense into R enviroment. Double_t RosenBrock(const TVectorD xx ) { const Double_t x = xx[0]; const Double_t y = xx[1]; const Double_t tmp1 = y-x*x; const Double_t tmp2 = 1-x; return 100*tmp1*tmp1+tmp2*tmp2; } TVectorD RosenBrockGrad(const TVectorD xx ) { const Double_t x = xx[0]; const Double_t y = xx[1]; TVectorD grad(2); grad[0]=-400 * x * (y - x * x) - 2 * (1 - x); grad[1]=200 * (y - x * x); return grad; } void Minimization() { #if defined(__CINT__) && !defined(__MAKECINT__) cout << "WARNING: This tutorial can run only using ACliC, you must run it by doing: " << endl; cout << "cd $ROOTSYS/tutorials/r/" << endl; cout << "\t .x Minimization.C+" << endl; return; #endif ROOT::R::TRInterface &r=gR->Instance(); //passsing RosenBrock funtion to R r["RosenBrock"]=ROOT::R::TRFunction(RosenBrock); //passsing RosenBrockGrad funtion to R r["RosenBrockGrad"]=ROOT::R::TRFunction(RosenBrockGrad); //the option "method" could be "Nelder-Mead", "BFGS", "CG", "L-BFGS-B", "SANN","Brent" //the option "control" let you put some constraints like //"maxit" The maximum number of iterations //"abstol" The absolute convergence tolerance. //"reltol" Relative convergence tolerance. r.Parse("result <- optim( c(0.01,0.01), RosenBrock,method='BFGS',control = list(maxit = 1000000) )"); //Getting result from R TVectorD min=r.ParseEval("result$par").ToVector<Double_t>(); std::cout.precision(8); //printing results std::cout<<"-----------------------------------------"<<std::endl; std::cout<<"Minimum x="<<min[0]<<" y="<<min[1]<<std::endl; std::cout<<"Value at minimum ="<<RosenBrock(min)<<std::endl; //using the gradient r.Parse("optimHess(result$par, RosenBrock, RosenBrockGrad)"); r.Parse("hresult <- optim(c(-1.2,1), RosenBrock, NULL, method = 'BFGS', hessian = TRUE)"); //getting the min calculated with the gradient TVectorD hmin=r.ParseEval("hresult$par").ToVector<Double_t>(); //printing results std::cout<<"-----------------------------------------"<<std::endl; std::cout<<"Minimization with the Gradient"<<endl; std::cout<<"Minimum x="<<hmin[0]<<" y="<<hmin[1]<<std::endl; std::cout<<"Value at minimum ="<<RosenBrock(hmin)<<std::endl; }
#include<TRInterface.h> #include<TRandom.h> void Interpolation() { //creting points TRandom r; TVectorD x(10),y(10); for(int i=0;i<10;i++) { x[i]=i; y[i]=r.Gaus(); } (*gR)["x"]=x; (*gR)["y"]=y; gR->Xwin();//required to active new window for plot //Plot parameter. Plotting using two rows in one column gR->Parse("par(mfrow = c(2,1))"); //plotting the points gR->plot("x, y, main = 'approx(.) and approxfun(.)'"); //approx returns a list with components x and y, //containing n coordinates which interpolate the given data points according to the method (and rule) desired. gR->points("approx(x, y), col = 2, pch = '*'"); gR->points("approx(x, y, method = 'constant'), col = 4, pch = '*'"); //The function approxfun returns a function performing (linear or constant) //interpolation of the given data points. //For a given set of x values, this function will return the corresponding interpolated values. gR->Parse("f <- approxfun(x, y)"); gR->curve("f(x), 0, 11, col = 'green2'"); gR->points("x, y"); //using approxfun with const method gR->Parse("fc <- approxfun(x, y, method = 'const')"); gR->curve("fc(x), 0, 10, col = 'darkblue', add = TRUE"); // different extrapolation on left and right side : gR->plot("approxfun(x, y, rule = 2:1), 0, 11,col = 'tomato', add = TRUE, lty = 3, lwd = 2"); }
#include<TRInterface.h> gR->Install("DEoptim");
//Example based in //http://cran.r-project.org/web/packages/DEoptim/DEoptim.pdf //Author: Omar Zapata #include<TRInterface.h> #include<TBenchmark.h> #include<math.h> #include<stdlib.h> //in the next function the pointer *double should be change by TVectorD, because the pointer has no //sense into R enviroment. It is a Generalization or RosenBrock Double_t GenRosenBrock(const TVectorD xx ) { int length=xx.GetNoElements(); Double_t result=0; for(int i=0;i<(length-1);i++) { result+=pow(1-xx[i],2)+100*pow(xx[i+1]-pow(xx[i],2),2); } return result; } void GlobalMinimization() { #if defined(__CINT__) && !defined(__MAKECINT__) cout << "WARNING: This tutorial can run only using ACliC, you must run it by doing: " << endl; cout << "cd $ROOTSYS/tutorials/r/" << endl; cout << "\t .x GlobalMinimization.C+" << endl; return; #endif TBenchmark bench; ROOT::R::TRInterface &r=gR->Instance(); //loading DEoptim r.Parse("require(DEoptim)"); // passsing RosenBrock funtion to R r["GenRosenBrock"]=ROOT::R::TRFunction(GenRosenBrock); //maximun of iterations r["MaxIter"]=5000; //n = of argument(vector) for GenRosenBrock r["n"]=10; //lower limits r.Parse("ll<-rep(-25, n)"); //upper limits r.Parse("ul<-rep(25, n)"); bench.Start("GlobalMinimization"); //calling minimization in sequential with one core and timing it. r.Parse("result<-DEoptim(fn=GenRosenBrock,lower=ll,upper=ul,control=list(NP=10*n,itermax=MaxIter))"); r.Parse("result"); std::cout<<"-----------------------------------------"<<std::endl; std::cout<<"Bechmark Times"<<std::endl; // printing times bench.Show("GlobalMinimization"); }
El documento original está disponible en https://clustercien.udea.edu.co/web/tiki-index.php?page=ROOTR+Cint