next up previous contents
Next: About this document Up: Beginner's Guide to MPI Previous: Sample Makefile

MPI Program with Graphics - Mandelbrot Rendering

 

This program generates a rendering of the Mandelbrot set in parallel.

What's the Mandelbrot set?

The Mandelbrot set is what is known as an iterative fractal. The idea is that every point (x, y) on the plane is mapped to a complex number c = x + yi. Using this value of c, we define a function f:

f(z) = z*z + c

We then define a series in which the nth element is f composed with itself n - 1 times applied to 0.

e.g., 0, f(0), f(f(0)), f(f(f(0))), f(f(f(f(0)))), ...

If this series is bounded, then point c is in the Mandelbrot set. If the series is not bounded, then c is not in the Mandelbrot set. It turns out that if the distance from any element in the series to the origin is greater than 2.0, then the series is not bounded. Thus, we can find out if a point is in the Mandelbrot set by iterating until the generated value is more than 2.0 away from the origin. Of course, if the point is in the set, then the iteration will never stop. Thus, we define an iteration limit N. Any point which can be iterated N times without escaping the radius 2.0 is probably in the Mandelbrot set.

Traditionally, images of the Mandelbrot set are drawn with the set in black. Areas outside the set are colored according to how many iterations there were before the value escaped radius 2.0. Such images are of course very pretty, and way overused.

Parallelizing Mandelbrot Rendering

How should you divide up the task of rendering the Mandelbrot set for N processors? The obvious method is to divide the image up into N sections, and have each processor render one section. This is inefficient, though, since different areas of the image take different amounts of time to compute, so some nodes will be done before others and sit idle. Even if the sections did require equal amounts of computation, some nodes may be faster than others.

The solution is to use manager-worker techniques. The idea is that one node (presumably node 0) is the manager, and farms out tasks to the rest of the nodes, which are the workers. The manager keeps a list of the workers, and whether each worker is busy. Whenever the manager needs a task done, it sends a message to an idle worker. When the worker completes, it sends a message back to the manager saying it is finished. Thus the manager can keep all the workers busy most of the time.

User Interface

When the program is started, it should bring up a graphics window containing a view of the entire Mandelbrot set. The user should be able to use the mouse to zoom in on a region. To quit the program, click the right mouse button inside the graphics window.

globals.h

#ifndef __globals_h__
#define __globals_h__

/* Includes: */
#include <stdlib.h>
#include <stdio.h>
#include <mpi.h>
#include "/usr/local/mpi/mpe/mpe.h"

/* Constants: */
enum {
  defaultWidth = 400,
  defaultHeight = 400
};

/* Variables: */
int rank, numnodes;
MPE_XGraph w;
int width, height;
int numcolors;
MPE_Color *colors;

/* Types: */
typedef double coordinate;

typedef struct {
  coordinate x1, y1, x2, y2;
  int left, top, right, bottom;
  int maxiter;
} MandelRender;

extern MPI_Datatype MandelRenderType;

/* Message tags: */
enum {
  tag_renderthis,
  tag_donerendering,
  tag_flush,
  tag_shutdown
};

#endif

mandel.c

#include "globals.h"
#include "manager.h"
#include "worker.h"

MPI_Datatype MandelRenderType;

void Setup_Types (void)
{
  int blockcounts[] = { 4, 5 };
  MPI_Aint offsets[] = { 0, 4 * sizeof(double) };
  MPI_Datatype types[2];
  types[0] = MPI_DOUBLE;
  types[1] = MPI_INT;
  MPI_Type_struct (2, blockcounts, offsets, types, &MandelRenderType);
  MPI_Type_commit (&MandelRenderType);
}

int main (int argc, char **argv)
{
  MPI_Init (&argc, &argv);
  MPI_Comm_rank (MPI_COMM_WORLD, &rank);
  MPI_Comm_size (MPI_COMM_WORLD, &numnodes);
  Setup_Types ();

  /* Open the window: */
  width = defaultWidth;
  height = defaultHeight;
  MPE_Open_graphics (&w, MPI_COMM_WORLD, NULL, -1, -1, 
		     width, height, 0);
  numcolors = 32;
  colors = (MPE_Color *) malloc(numcolors * sizeof(MPE_Color));
  MPE_Make_color_array (w, numcolors, colors);
  
  /* Start the manager or worker, as appropriate: */
  if (!rank)
    Manager ();
  else
    Worker ();

  MPE_Close_graphics (&w);

  MPI_Finalize ();
  exit(0);
}

manager.h

#include "globals.h"
void Manager (void);

manager.c

#include "/usr/local/users/hisley/mpi/src/lab4/manager.h"
#include "/usr/local/users/hisley/mpi/src/lab4/mouse_status.h"

/* The busy array contains an element for every worker node.  
   The value of node[n] indicates how many tasks node n is doing.
   We try to keep each worker fed with 2 tasks at a time, to avoid message 
   latency.  Thus workersfree specifies how many more tasks could be fed to
   workers. */
char *workersbusy;
int workersfree;

static int ReceiveDonePacket (void)
{
  MPI_Status status;
  MPI_Recv (NULL, 0, MPI_INT, MPI_ANY_SOURCE, tag_donerendering, 
	    MPI_COMM_WORLD, &status);
  workersbusy[status.MPI_SOURCE]--;
  workersfree++;
  return status.MPI_SOURCE;
}

static void Manager_Draw (MandelRender *r)
{
  int sectHeight;
  coordinate sectDeltaY;
  MandelRender section;
  int who;

  /* Calculate height of the sections to render: */
  sectHeight = height / (numnodes - 1) / 10;	/* Let's go for 10 times as 
						   many sections as workers. */
  if (sectHeight <= 0)
    sectHeight = 1;
  sectDeltaY = (r->y2 - r->y1) / (r->bottom - r->top) * sectHeight;
  
  section = *r;
  while (1) {
    section.bottom = section.top + sectHeight;
    section.y2 = section.y1 + sectDeltaY;
    if (section.bottom > r->bottom) {
      section.bottom = r->bottom;
      section.y2 = r->y2;
    }
    
    if (workersfree)
      for (who = 1; who < numnodes && workersbusy[who] == 2; who++);
    else
      who = ReceiveDonePacket ();
      
    MPI_Send (&section, 1, MandelRenderType, who, tag_renderthis,
	      MPI_COMM_WORLD);
    workersbusy[who]++;
    workersfree--;

    if (section.bottom == r->bottom)
      break;
    
    section.top = section.bottom;
    section.y1 = section.y2;
  }

  for (who = 1; who < numnodes; who++)
    MPI_Send (NULL, 0, MPI_INT, who, tag_flush, MPI_COMM_WORLD);

  while (workersfree != (numnodes - 1) * 2)
    ReceiveDonePacket ();
}

void Manager (void)
{
  MandelRender whole = { -2.0, 1.5, 1.0, -1.5,
			 0, 0, -1, -1, 1000 },
	       r;
  char done = 0;

  puts ("Welcome to the lame Mandelbrot rendering demo.\n\n"
	"When I say I'm ready for mouse input you may:\n"
	"  - Drag out a square region to zoom in on with the left mouse "
	"button.\n"
	"  - Click the middle button to return to a view of the full "
	"Mandelbrot set.\n"
	"  - Click the right button to quit this program.\n");

  whole.right = width;
  whole.bottom = height;
  r = whole;
  workersbusy = (char *) calloc (sizeof(char), numnodes);
  workersfree = (numnodes - 1) * 2;

  while (!done) {
    printf ("\nDrawing image:\nReal: %12.10f - %12.10f; 
             Imaginary: %12.10f - %12.10f.\n", r.x1, r.x2, r.y1, r.y2);
    fflush (stdout);
    Manager_Draw (&r);
    while (1) {
      int button, x, y;
      printf("\nReady for mouse input.\n");
      fflush (stdout);
      MPE_Get_mouse_press (w, &x, &y, &button);
      if (button == 3) {
	done = 1;
	break;
      } else if (button == 2) {
	r = whole;
	break;
      } else if (button == 1) {
	int rx, ry;
	double ox1 = r.x1, oy1 = r.y1;
	MPE_Drag_square (w, &x, &y, &rx, &ry); 
	r.x1 = ox1 + (r.x2 - ox1) * x / width;
	r.x2 = ox1 + (r.x2 - ox1) * rx / width;
	r.y1 = oy1 + (r.y2 - oy1) * y / height;
	r.y2 = oy1 + (r.y2 - oy1) * ry / height;
	break;
      }
    }
  }

  { int who;
    for (who = 1; who < numnodes; who++)
      MPI_Send (NULL, 0, MPI_INT, who, tag_shutdown, MPI_COMM_WORLD);
  }
}

worker.h

#include "/usr/local/users/hisley/mpi/src/lab4/globals.h"
void Worker (void);

worker.c

#include "/usr/local/users/hisley/mpi/src/lab4/worker.h"

static void RenderSection (MandelRender *r)
{
  int width = r->right - r->left, height = r->bottom - r->top;
  int numpoints = width * height;
  int x, y;
  coordinate dr = (r->x2 - r->x1) / (r->right - r->left),
	     di = (r->y2 - r->y1) / (r->bottom - r->top),
	     cr, ci = r->y1;
  MPE_Point *points = (MPE_Point *) malloc (numpoints * sizeof(MPE_Point)),
	    *point = points;
  
  for (y = r->top; y < r->bottom; y++, (ci += di)) {
    cr = r->x1;
    for (x = r->left; x < r->right; x++, (cr += dr)) {
      int iter;
      register double zr = 0, zi = 0;
      point->x = x;
      point->y = y;
      for (iter = r->maxiter; iter && zr * zr + zi * zi < 4.0; iter--) {
	register double nr = zr * zr - zi * zi + cr;
	zi = zr * zi;
	zi += zi + ci;
	zr = nr;
      }
      point->c = iter ? colors[iter % numcolors] : MPE_BLACK;
      point++;
    }
  }

  MPE_Draw_points (w, points, numpoints);
  free (points);
}

void Worker (void)
{
  MandelRender r;
  MPI_Status status;
  while (1) {
    MPI_Probe (0, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
    if (status.MPI_TAG == tag_renderthis) {
      MPI_Recv (&r, 1, MandelRenderType, status.MPI_SOURCE, tag_renderthis, 
		MPI_COMM_WORLD, &status);
      RenderSection (&r);
      MPI_Send (NULL, 0, MPI_INT, status.MPI_SOURCE, tag_donerendering,
		MPI_COMM_WORLD);
    } else if (status.MPI_TAG == tag_flush) {
      MPI_Recv (NULL, 0, MPI_INT, status.MPI_SOURCE, tag_flush, 
               MPI_COMM_WORLD, &status);
      MPE_Update (w);
    } else if (status.MPI_TAG == tag_shutdown) {
      MPI_Recv (NULL, 0, MPI_INT, status.MPI_SOURCE, tag_shutdown, 
               MPI_COMM_WORLD, &status);
      return;
    } else {
      fprintf (stderr, "Hey!  Unknown tag %d on node %d.\n", 
	      status.MPI_TAG, rank);
      fflush (stderr);
    }
  }
}

mouse_status.h

#ifndef __mouse_status_h__
#define __mouse_status_h__

#include "/usr/local/mpi/mpe/mpe.h"

int MPE_Get_mouse_status( MPE_XGraph graph, int *x, int *y, 
			  int *button, int *wasPressed );

void MPE_Drag_square (MPE_XGraph w, int *startx, int *starty, 
		      int *endx, int *endy);

#endif

mouse_status.c

#include <math.h>
#include "/usr/local/mpi/mpe/basex11.h"

#define MPE_INTERNAL
#include "/usr/local/mpi/mpe/mpe.h"   /*I "mpe.h" I*/

#include "/usr/local/users/hisley/mpi/src/lab4/mouse_status.h"

#define DEBUG 0

/*@
  MPE_Get_mouse_status - Checks for mouse button location & stuff
  Checks if the mouse button has been pressed inside this MPE window.
  If pressed, returns the coordinate relative to the upper right of
  this MPE window and the button that was pressed.

  Input Parameter:
. graph - MPE graphics handle

  Output Parameters:
. x - horizontal coordinate of the point where the mouse button was pressed
. y - vertical coordinate of the point where the mouse button was pressed
. button - which button was pressed: MPE_BUTTON[1-5]
. wasPressed - 1 if the button was pressed, 0 if not

@*/
int MPE_Get_mouse_status( MPE_XGraph graph, int *x, int *y, 
			  int *button, int *wasPressed )
{
  Window root, child;
  int rx, ry;
  unsigned int mask;

  if (graph->Cookie != MPE_G_COOKIE) {
    fprintf( stderr, "Handle argument is incorrect or corrupted\n" );
    return MPE_ERR_BAD_ARGS;
  }

  XQueryPointer (graph->xwin->disp, graph->xwin->win, &root, &child, &rx, &ry,
		 x, y, &mask);
  if (mask & Button1Mask)
    *button = Button1;
  else if (mask & Button2Mask)
    *button = Button2;
  else if (mask & Button3Mask)
    *button = Button3;
  else
    *button = 0;
  *wasPressed = *button;

  return MPE_SUCCESS;
}


static void fill_ones_rect (MPE_XGraph handle, int x, int y, int w, int h)
{
  XBSetPixVal (handle->xwin, 0xFFFFFFFF);
  XFillRectangle (handle->xwin->disp, handle->xwin->win,
		  handle->xwin->gc.set, x, y, w, h);
}


void MPE_Drag_square (MPE_XGraph w, int *startx, int *starty, 
		      int *endx, int *endy)
{
  int button, pressed, x, y, ox = *startx, oy = *starty;

  MPE_Draw_logic (w, MPE_LOGIC_XOR);
  MPE_Get_mouse_status (w, &x, &y, &button, &pressed);
  if (!pressed)
    MPE_Get_mouse_press (w, &x, &y, &button);
  while (pressed) {
    MPE_Get_mouse_status (w, &x, &y, &button, &pressed);
    y = *starty + (x - *startx);
    if (x != ox) {
      if (ox > *startx)
	fill_ones_rect (w, *startx, *starty, ox - *startx, oy - *starty);
      else
	fill_ones_rect (w, ox, oy, *startx - ox, *starty - oy);
      ox = x; oy = y;
      if (pressed)
	if (ox > *startx)
	  fill_ones_rect (w, *startx, *starty, ox - *startx, oy - *starty);
	else
	  fill_ones_rect (w, ox, oy, *startx - ox, *starty - oy);
      MPE_Update (w);
    }
  }
  if (x > *startx) {
    *endx = x; *endy = y;
  } else {
    *endx = *startx; *endy = *starty;
    *startx = x; *starty = y;
  }
    
  MPE_Draw_logic (w, MPE_LOGIC_COPY);
}



next up previous contents
Next: About this document Up: Beginner's Guide to MPI Previous: Sample Makefile



Lori Pollock
Wed Feb 4 14:18:58 EST 1998