Cancel Fullscreen
Imprimir Comparte esta página

RootMpi Demos

Image


Index




Go to Index

Sphere


In this demonstration I pretend show you how to think a parallel problem and the power of to use RootMpi to solve it using 4 processors exactly.
I want to calculate the parallel surface and volume of a sphere.

Parallel Sphere Volume
To calculate Volume I can calculate $ V=\int\limits_{0}^{r} r^2-(x^2)\,dx. $ . I can calculate it using ROOT::Mpi::Math::TIntegratorWrapper that is RootMpi's class with parallelization ready.
/*********************************************************************/
   /*Section to calculate volume with parallel multimensional integrator*/
   /*********************************************************************/   
   //integrator creation root's style
   ROOT::Math::Functor1D function(&sphere_function);

   ROOT::Math::Integrator integrator;
   integrator.SetFunction(function);

   //Mpi's integrator wrapper, integrator from ROOT in one dimension
   ROOT::Mpi::Math::TIntegratorWrapper MpiIntegrator(&integrator);

   //integrate in parallel and send the result to rank zero, intervar [0,radius]
   MpiIntegrator.Integrate(0, 0, radius);


Parallel Sphere Surface
(wait 5 secs to render latex)
To calculate the surface I can express spherical coordinates in rectangular coordinates in this way:
  • $x = rcos(\theta)sin(\phi)$
  • $y = rsin(\theta)sin(\phi)$
  • $z = rcos(\phi)$
$\phi \in [0,2\pi] and\ \theta \in [0,\pi]$
Range for process
Process 0 $\phi \in [0,0.5\pi] and\ \theta \in [0,\pi]$
Process 1 $\phi \in [0.5,\pi] and\ \theta \in [0,\pi]$
Process 2 $\phi \in [1,1.5\pi] and\ \theta \in [0,\pi]$
Process 3 $\phi \in [1.5,2\pi] and\ \theta \in [0,\pi]$
In every section I evaluate the intervals with delta=0.01, then I have almost 50000 points in every process.
Every color is the region calculate in different processor, the points are encapsulated in a object TPolyMarker3D and sent to root process(zero) and painted.
Reference(external link)


Image


This is the macro for ROOT interpreter.

Go to Index
SphereCode

Copy To Clipboard
  1 . /*******************************************************************************/
  2 . /*Author: Omar Andres Zapata Mesa  Ago 2012                                    */
  3 . /*this example calculate parallel surface's regions in a sphere and the volume.*/
  4 . /*need 4 process to work, beacase the sphere was splitted into four parts.     */
  5 . /*******************************************************************************/
  6 . #include <Rtypes.h>
  7 . #include <Math/Math.h>
  8 . #include <Math/Functor.h>
  9 . #include <Math/Integrator.h>
 10 . 
 11 . gSystem->Load("libRMpi");
 12 . 
 13 . Double_t delta=0.01;//delta to iterate points in the surface
 14 . Double_t radius=10; //sphere radius
 15 . Double_t radius2=radius*radius; //radius' square
 16 . Double_t theta, phi, x,y,z;
 17 . Int_t point=0;
 18 . /************************************************************************************/
 19 . /*this is the function for positive region of a sphere. use to integrate the volume.*/
 20 . /************************************************************************************/
 21 . Double_t sphere_function(const Double_t param)
 22 . {
 23 .   Double_t x,y;
 24 .   x=param;//coor x
 25 .   return 2*M_PI*(radius2-(x*x));//funtion sphere_function(r,x)=V(r,x)= 2*pi*(r2-x2)
 26 . }
 27 . 
 28 . void sphere()
 29 . {
 30 .    /********************/
 31 .    /*Mpi initialization*/
 32 .    /********************/
 33 .    ROOT::Mpi::TEnvironment env(gApplication->Argc(), gApplication->Argv());
 34 .    ROOT::Mpi::TIntraCommunicator world;
 35 . 
 36 .    /*********************************/
 37 .    /*Error control, need 4 processes*/
 38 .    /*********************************/   
 39 .    if(world.Size()!=4)
 40 .    {
 41 .       std::cout<<"Error: You need run with 4 processes exactly, the sphere was spplitted in four regions. \nExample: mpirun -np 4 root -l  sphere.C\n";
 42 .       std::cout.flush();
 43 .       world.Abort(6);  
 44 .    }
 45 .    
 46 .    /*********************************************************************/
 47 .    /*Section to calculate volume with parallel multimensional integrator*/
 48 .    /*********************************************************************/   
 49 .    //integrator creation root's style
 50 .    ROOT::Math::Functor1D function(&sphere_function);
 51 . 
 52 .    ROOT::Math::Integrator integrator;
 53 .    integrator.SetFunction(function);
 54 . 
 55 .    //Mpi's integrator wrapper, integrator from ROOT in one dimension
 56 .    ROOT::Mpi::Math::TIntegratorWrapper MpiIntegrator(&integrator);
 57 . 
 58 .    //integrate in parallel and send the result to rank zero, intervar [0,radius]
 59 .    MpiIntegrator.Integrate(0, 0, radius);
 60 . 
 61 .    
 62 .    /*****************************************************************************************/
 63 .    /*this is the root process,calculate the surface z=+qsrt(1−(x2+b2)) for x∈[−1,0],y∈[−1,1]*/
 64 .    /*and recv the other surface regions in TPolyMarker3D object from others ranks           */
 65 .    /*****************************************************************************************/
 66 .    if (world.Rank() == 0) {
 67 .       gBenchmark->Start("sphere");
 68 .       // create and open a canvas
 69 .       TCanvas *canvas = new TCanvas("sphere", "sphere", 300, 10, 700, 500);
 70 .       canvas->SetFillColor(14);
 71 .       // creating view
 72 .       TView *view = TView::CreateView(1, 0, 0);
 73 .       view->RotateView(90,-10);
 74 .       view->SetRange(-radius*2, -radius*2, -2*radius, 2*radius, 2*radius, 2*radius);
 75 .       // create a PolyMarker3D for region/segment 1 calculate in rank 0(this rank)
 76 .       TPolyMarker3D *surfaceSeg13d = new TPolyMarker3D();
 77 .       
 78 .       for (theta = 0; theta <0.5*M_PI ; theta+= delta) {
 79 .          for (phi = 0; phi <M_PI; phi += delta) {
 80 . 	   
 81 . 	   x=radius*cos(theta)*sin(phi);
 82 . 	   y=radius*sin(theta)*sin(phi);
 83 . 	   z=radius*cos(phi);
 84 .             // set points
 85 .             surfaceSeg13d->SetPoint(point, x, y, z);
 86 .             point++;
 87 .          }
 88 .       }
 89 .       // set marker size, color & style
 90 .       surfaceSeg13d->SetMarkerSize(1);
 91 .       surfaceSeg13d->SetMarkerColor(world.Rank() + world.Size());
 92 .       surfaceSeg13d->SetMarkerStyle(1);
 93 .       //drawing region 1
 94 .       surfaceSeg13d->Draw();
 95 . 
 96 .       //getting region 2 from rank 1
 97 .       TPolyMarker3D *surfaceSeg23d = new TPolyMarker3D;
 98 .       world.RecvObject(1, 1, *surfaceSeg23d);
 99 .       surfaceSeg23d->Draw();
100 . 
101 .       //getting region 3 from rank 2
102 .       TPolyMarker3D *surfaceSeg33d = new TPolyMarker3D;
103 .       world.RecvObject(2, 1, *surfaceSeg33d);
104 .       surfaceSeg33d->Draw();
105 . 
106 .       //getting region 4 from rank 3
107 .       TPolyMarker3D *surfaceSeg43d = new TPolyMarker3D;
108 .       world.RecvObject(3, 1, *surfaceSeg43d);
109 .       surfaceSeg43d->Draw();
110 . 
111 .       Char_t timeStr[60];
112 .       Char_t pointsStr[100];
113 .       Char_t radiusStr[100];
114 .       Char_t deltaStr[100];
115 .       gBenchmark->Show("sphere");
116 . 
117 .       Float_t ct = gBenchmark->GetCpuTime("sphere");
118 .       sprintf(timeStr, "Execution time: %g sec.", ct);
119 .       sprintf(radiusStr, "Radius : %g", radius);
120 .       sprintf(deltaStr, "Delta : %g", delta);
121 .       sprintf(pointsStr, "Total surface points : %g", point*4);
122 .       
123 .       TPaveText *text = new TPaveText(0.1, 0.6, 0.9, 0.99);
124 .       text->SetFillColor(42);
125 .       text->AddText("ROOTMpi example: sphere.C");
126 .       text->AddText(timeStr);
127 .       text->AddText(radiusStr);
128 .       text->AddText(deltaStr);
129 .       text->AddText(pointsStr);
130 .       text->Draw();
131 .       
132 .       //sphere volume  V=(4/3)*pi*r3
133 .       //and compare the result integration's volume
134 .       Char_t volumeStr[100];
135 .       Double_t volume= MpiIntegrator.Result();
136 .       Double_t volumeTeo=(4/3.0)*M_PI*(radius*radius*radius);
137 .       
138 .       sprintf(volumeStr, "Volumen with integration = %.15g \nVolumen with (4/3)*pi*r3 = %.15g", volume,volumeTeo);
139 .       
140 .       TPaveText *volumeText = new TPaveText(-0.1, -0.81, -0.9, -0.97);
141 .       volumeText->SetFillColor(41);
142 .       volumeText->SetTextAlign(11);
143 .       volumeText->AddText(volumeStr);
144 .       volumeText->Draw();
145 .    }
146 .    /*************************************************************************/
147 .    /*this rank calculate the surface for theta ∈[0.5pi,pi] and phi ∈[0,pi]  */
148 .    /*and send the surface regions in TPolyMarker3D object to rank root zero */
149 .    /*************************************************************************/
150 .    if (world.Rank() == 1) {
151 .       TPolyMarker3D *surfaceSeg23d = new TPolyMarker3D();
152 . 
153 .       for (theta = 0.5*M_PI; theta <M_PI ; theta+= delta) {
154 .          for (phi = 0; phi <M_PI; phi += delta) {
155 . 	   
156 . 	   x=radius*cos(theta)*sin(phi);
157 . 	   y=radius*sin(theta)*sin(phi);
158 . 	   z=radius*cos(phi);
159 .             // set points
160 .             surfaceSeg23d->SetPoint(point, x, y, z);
161 .             point++;
162 .          }
163 .       }
164 .       // set marker size, color & style
165 .       surfaceSeg23d->SetMarkerSize(1);
166 .       surfaceSeg23d->SetMarkerColor(world.Rank() + world.Size());
167 .       surfaceSeg23d->SetMarkerStyle(1);
168 .       //send surface points to root process
169 .       world.SendObject(0, 1, *surfaceSeg23d);
170 .    }
171 .    /*************************************************************************/
172 .    /*this rank calculate the surface for theta ∈[pi,1.5pi] and phi ∈[0,pi]  */
173 .    /*and send the surface regions in TPolyMarker3D object to rank root zero */
174 .    /*************************************************************************/
175 .    if (world.Rank() == 2) {
176 .       TPolyMarker3D *surfaceSeg33d = new TPolyMarker3D();
177 .       for (theta = M_PI; theta <1.5*M_PI ; theta+= delta) {
178 .          for (phi = 0; phi <M_PI; phi += delta) {
179 . 	   
180 . 	   x=radius*cos(theta)*sin(phi);
181 . 	   y=radius*sin(theta)*sin(phi);
182 . 	   z=radius*cos(phi);
183 .             // set points
184 .             surfaceSeg33d->SetPoint(point, x, y, z);
185 .             point++;
186 .          }
187 .       }
188 .       // set marker size, color & style
189 .       surfaceSeg33d->SetMarkerSize(1);
190 .       surfaceSeg33d->SetMarkerColor(world.Rank() + world.Size());
191 .       surfaceSeg33d->SetMarkerStyle(1);
192 .       //send surface points to root process
193 .       world.SendObject(0, 1, *surfaceSeg33d);
194 .    }
195 .    /*************************************************************************/
196 .    /*this rank calculate the surface for theta ∈[1.5pi,2pi] and phi ∈[0,pi]  */
197 .    /*and send the surface regions in TPolyMarker3D object to rank root zero */
198 .    /*************************************************************************/
199 .    if (world.Rank() == 3) {
200 .       TPolyMarker3D *surfaceSeg43d = new TPolyMarker3D();
201 .       for (theta = 1.5*M_PI; theta <2*M_PI ; theta+= delta) {
202 .          for (phi = 0; phi <M_PI; phi += delta) {
203 . 	   
204 . 	   x=radius*cos(theta)*sin(phi);
205 . 	   y=radius*sin(theta)*sin(phi);
206 . 	   z=radius*cos(phi);
207 .             // set points
208 .             surfaceSeg43d->SetPoint(point, x, y, z);
209 .             point++;
210 .          }
211 .       }
212 .       // set marker size, color & style
213 .       surfaceSeg43d->SetMarkerSize(1);
214 .       surfaceSeg43d->SetMarkerColor(world.Rank() + world.Size());
215 .       surfaceSeg43d->SetMarkerStyle(1);
216 .       //send surface points to root process
217 .       world.SendObject(0, 1, *surfaceSeg43d);
218 .    }
219 . }


Go to Index
SphereVideo

Created by Omar Andres Zapata Mesa.
Voice by Gerardo Gutierrez. (Thanks!)


Go to Index

Parallel Wave Animation

In this demonstration we are going to animate a sin wave. We split it into n-1 processes, the process root(zero) gets the points from other processes and it sends a broadcast with the new angle to calculate the height, and repeat the procedure.
WARNING this animation require OpenGL support in your pc.


Image

Go to Index
WaveCode

  1 . /*******************************************************************************/
  2 . /*Author: Omar Andres Zapata Mesa  Ago 2012                                    */
  3 . /*this example calculate parallel wave animation                               */
  4 . /*work with any number(n) of processes but the wave is splitted into n-1 parts */
  5 . /*******************************************************************************/
  6 . #include <Rtypes.h>
  7 . #include <Math/Math.h>
  8 . #include <TF1.h>
  9 . #include <TPolyMarker3D.h>
 10 . Double_t delta=0.1;
 11 . Double_t min=0,max=6*M_PI;
 12 . Double_t angle=M_PI;
 13 . Double_t init,end;
 14 . Int_t counter=0;
 15 . 
 16 . UInt_t rep=1000;
 17 . UInt_t delay=100;
 18 . 
 19 . UInt_t MPI_SUCCESS=0;
 20 . /********************/
 21 . /*Mpi initialization*/
 22 . /********************/
 23 . ROOT::Mpi::TEnvironment env(gApplication->Argc(), gApplication->Argv());
 24 . ROOT::Mpi::TIntraCommunicator world;
 25 . 
 26 . TCanvas *canvas = NULL;
 27 . TView   *view   = NULL;
 28 . 
 29 . 
 30 . void Init()
 31 . {
 32 .  gStyle->SetCanvasPreferGL(true);
 33 .  gStyle->SetFrameFillColor(42);
 34 .  canvas = new TCanvas("canvas", "wave", 300, 10, 700, 500);
 35 .  canvas->SetFillColor(14);
 36 .  view = TView::CreateView(1, 0, 0);
 37 .  view->RotateView(90,0);
 38 . }
 39 . 
 40 . void Animate()
 41 . {
 42 .      if (!gROOT->GetListOfCanvases()->FindObject("canvas")){
 43 .        world.Abort(MPI_SUCCESS);
 44 .     };
 45 .      gPad->Clear();
 46 .      for(int i=1;i<world.Size();i++)
 47 .       {
 48 .       //getting wave sections from ranks>1
 49 .       TPolyMarker3D *Section = new TPolyMarker3D();
 50 .       world.RecvObject(i, 1, *Section);
 51 .       Section->Draw();
 52 .       }
 53 .       gPad->Modified();
 54 .       gPad->Update();
 55 .       angle += 0.05*M_PI;
 56 .       world.BcastSendScalar(angle,0);
 57 .       counter++;
 58 .       if(counter<rep) world.BcastSendScalar(true,0);
 59 .       else world.BcastSendScalar(false,0);
 60 . }
 61 . 
 62 . 
 63 . void wave()
 64 . {
 65 .    
 66 .    /**************************************************/
 67 .    /*rank zero(root) dont calculate wave, init viewer*/
 68 .    /*and recv wave sections from other processes     */
 69 .    /**************************************************/
 70 .    if (world.Rank() == 0) {
 71 .       Init();
 72 .       TTimer *timer = new TTimer(delay);
 73 .       timer->SetCommand("Animate()");
 74 .       Char_t timeStr[60];
 75 .       Char_t pointsStr[100];
 76 .       Char_t deltaStr[100];      
 77 .       timer->TurnOn();
 78 .   }else{
 79 .    Double_t step = fabs(min - max) / (world.Size()-1); 
 80 . 
 81 .    init = min + (world.Rank()-1) * step;
 82 .    end  = init + step;
 83 .    cout<<"rank="<<world.Rank()<<" init="<<init<<" end="<<end<<"\n";
 84 .    Int_t flag=true;
 85 .    while(flag){
 86 .      TPolyMarker3D *Section=new TPolyMarker3D();
 87 .      counter=0;
 88 .           for(Double_t t=init;t<end;t+=delta)
 89 .           {
 90 .                 Section->SetPoint(counter,t,sin(t+angle),0);
 91 .                 counter++;
 92 . 	  }
 93 .       // set marker size, color & style
 94 .       Section->SetMarkerSize(1);
 95 .       Section->SetMarkerColor(world.Rank());
 96 .       Section->SetMarkerStyle(8);
 97 .       //send segment of points to root process
 98 .       world.SendObject(0, 1, *Section);
 99 .       
100 .       world.BcastRecvScalar(angle,0); //recv angle to wave evolution
101 .       world.BcastRecvScalar(flag,0);//recv flag to stop
102 .       delete Section;
103 .       }
104 .        world.Abort(MPI_SUCCESS);
105 .   }
106 . }




Go to Index
WaveVideo

Created by Omar Andres Zapata Mesa.
Voice by Gerardo Gutierrez. (Thanks!)


Go to Index

Visits