int main(int argc, char *argv[]) { complex_t coord_point, julia_constant; double x_max, x_min, y_max, y_min, x_resolution, y_resolution; double divergent_limit; char file_message[160]; char filename[100]; int icount, imax_iterations; int ipixels_across, ipixels_down; int i, j, k, julia, alternate_equation; int imin, imax, jmin, jmax; int work[5]; header_t *result = NULL; /* make an integer array of size [N x M] to hold answers. */ int *in_grid_array, *out_grid_array = NULL; int numprocs; int namelen; char processor_name[MPI_MAX_PROCESSOR_NAME]; int num_colors; color_t *colors = NULL; MPI_Status status; int listener; int save_image = 0; int optval; int num_children = DEFAULT_NUM_SLAVES; int master = 1; MPI_Comm parent, *child_comm = NULL; MPI_Request *child_request = NULL; int error_code, error; char error_str[MPI_MAX_ERROR_STRING]; int length; int index; int pid; MPI_Init(&argc, &argv); MPI_Comm_size(MPI_COMM_WORLD, &numprocs); MPI_Comm_rank(MPI_COMM_WORLD, &myid); MPI_Comm_get_parent(&parent); MPI_Get_processor_name(processor_name, &namelen); pid = getpid(); if (parent == MPI_COMM_NULL) { if (numprocs > 1) { printf("Error: only one process allowed for the master.\n"); PrintUsage(); error_code = MPI_Abort(MPI_COMM_WORLD, -1); exit(error_code); } printf("Welcome to the Mandelbrot/Julia set explorer.\n"); master = 1; /* Get inputs-- region to view (must be within x/ymin to x/ymax, make sure xmax>xmin and ymax>ymin) and resolution (number of pixels along an edge, N x M, i.e. 256x256) */ read_mand_args(argc, argv, &imax_iterations, &ipixels_across, &ipixels_down, &x_min, &x_max, &y_min, &y_max, &julia, &julia_constant.real, &julia_constant.imaginary, &divergent_limit, &alternate_equation, filename, &num_colors, &use_stdin, &save_image, &num_children); check_mand_params(&imax_iterations, &ipixels_across, &ipixels_down, &x_min, &x_max, &y_min, &y_max, &divergent_limit, &num_children); if (julia == 1) /* we're doing a julia figure */ check_julia_params(&julia_constant.real, &julia_constant.imaginary); /* spawn slaves */ child_comm = (MPI_Comm*)malloc(num_children * sizeof(MPI_Comm)); child_request = (MPI_Request*)malloc(num_children * sizeof(MPI_Request)); result = (header_t*)malloc(num_children * sizeof(header_t)); if (child_comm == NULL || child_request == NULL || result == NULL) { printf("Error: unable to allocate an array of %d communicators, requests and work objects for the slaves.\n", num_children); error_code = MPI_Abort(MPI_COMM_WORLD, -1); exit(error_code); } printf("Spawning %d slaves.\n", num_children); for (i=0; i<num_children; i++) { error = MPI_Comm_spawn(argv[0], MPI_ARGV_NULL, 1, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &child_comm[i], &error_code); if (error != MPI_SUCCESS) { error_str[0] = '\0'; length = MPI_MAX_ERROR_STRING; MPI_Error_string(error, error_str, &length); printf("Error: MPI_Comm_spawn failed: %s\n", error_str); error_code = MPI_Abort(MPI_COMM_WORLD, -1); exit(error_code); } } /* send out parameters */ for (i=0; i<num_children; i++) { MPI_Send(&num_colors, 1, MPI_INT, 0, 0, child_comm[i]); MPI_Send(&imax_iterations, 1, MPI_INT, 0, 0, child_comm[i]); MPI_Send(&ipixels_across, 1, MPI_INT, 0, 0, child_comm[i]); MPI_Send(&ipixels_down, 1, MPI_INT, 0, 0, child_comm[i]); MPI_Send(&divergent_limit, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&julia, 1, MPI_INT, 0, 0, child_comm[i]); MPI_Send(&julia_constant.real, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&julia_constant.imaginary, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&alternate_equation, 1, MPI_INT, 0, 0, child_comm[i]); } } else { master = 0; MPI_Recv(&num_colors, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&imax_iterations, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&ipixels_across, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&ipixels_down, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&divergent_limit, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&julia, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&julia_constant.real, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&julia_constant.imaginary, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&alternate_equation, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); } if (master) { colors = malloc((num_colors+1)* sizeof(color_t)); if (colors == NULL) { MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } Make_color_array(num_colors, colors); colors[num_colors] = 0; /* add one on the top to avoid edge errors */ } /* allocate memory */ if ( (in_grid_array = (int *)calloc(ipixels_across * ipixels_down, sizeof(int))) == NULL) { printf("Memory allocation failed for data array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } if (master) { int istep, jstep; int i1[400], i2[400], j1[400], j2[400]; int ii, jj; struct sockaddr_in addr; int len; char line[1024], *token; if ( (out_grid_array = (int *)calloc(ipixels_across * ipixels_down, sizeof(int))) == NULL) { printf("Memory allocation failed for data array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } srand(getpid()); if (!use_stdin) { addr.sin_family = AF_INET; addr.sin_addr.s_addr = INADDR_ANY; addr.sin_port = htons(DEFAULT_PORT); listener = socket(AF_INET, SOCK_STREAM, 0); if (listener == -1) { printf("unable to create a listener socket.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } if (bind(listener, &addr, sizeof(addr)) == -1) { addr.sin_port = 0; if (bind(listener, &addr, sizeof(addr)) == -1) { printf("unable to create a listener socket.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } } if (listen(listener, 1) == -1) { printf("unable to listen.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } len = sizeof(addr); getsockname(listener, &addr, &len); printf("%s listening on port %d\n", processor_name, ntohs(addr.sin_port)); fflush(stdout); sock = accept(listener, NULL, NULL); if (sock == -1) { printf("unable to accept a socket connection.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } printf("accepted connection from visualization program.\n"); fflush(stdout); #ifdef HAVE_WINDOWS_H optval = 1; setsockopt(sock, IPPROTO_TCP, TCP_NODELAY, (char *)&optval, sizeof(optval)); #endif printf("sending image size to visualizer.\n"); sock_write(sock, &ipixels_across, sizeof(int)); sock_write(sock, &ipixels_down, sizeof(int)); sock_write(sock, &num_colors, sizeof(int)); sock_write(sock, &imax_iterations, sizeof(int)); } for (;;) { /* get x_min, x_max, y_min, and y_max */ if (use_stdin) { printf("input xmin ymin xmax ymax max_iter, (0 0 0 0 0 to quit):\n");fflush(stdout); fgets(line, 1024, stdin); printf("read <%s> from stdin\n", line);fflush(stdout); token = strtok(line, " \n"); x_min = atof(token); token = strtok(NULL, " \n"); y_min = atof(token); token = strtok(NULL, " \n"); x_max = atof(token); token = strtok(NULL, " \n"); y_max = atof(token); token = strtok(NULL, " \n"); imax_iterations = atoi(token); /*sscanf(line, "%g %g %g %g", &x_min, &y_min, &x_max, &y_max);*/ /*scanf("%g %g %g %g", &x_min, &y_min, &x_max, &y_max);*/ } else { printf("reading xmin,ymin,xmax,ymax.\n");fflush(stdout); sock_read(sock, &x_min, sizeof(double)); sock_read(sock, &y_min, sizeof(double)); sock_read(sock, &x_max, sizeof(double)); sock_read(sock, &y_max, sizeof(double)); sock_read(sock, &imax_iterations, sizeof(int)); } printf("x0,y0 = (%f, %f) x1,y1 = (%f,%f) max_iter = %d\n", x_min, y_min, x_max, y_max, imax_iterations);fflush(stdout); /*printf("sending the limits: (%f,%f)(%f,%f)\n", x_min, y_min, x_max, y_max);fflush(stdout);*/ /* let everyone know the limits */ for (i=0; i<num_children; i++) { MPI_Send(&x_min, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&x_max, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&y_min, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&y_max, 1, MPI_DOUBLE, 0, 0, child_comm[i]); MPI_Send(&imax_iterations, 1, MPI_INT, 0, 0, child_comm[i]); } /* check for the end condition */ if (x_min == x_max && y_min == y_max) { /*printf("root bailing.\n");fflush(stdout);*/ break; } /* break the work up into 400 pieces */ istep = ipixels_across / 20; jstep = ipixels_down / 20; if (istep < 1) istep = 1; if (jstep < 1) jstep = 1; k = 0; for (i=0; i<20; i++) { for (j=0; j<20; j++) { i1[k] = MIN(istep * i, ipixels_across - 1); i2[k] = MIN((istep * (i+1)) - 1, ipixels_across - 1); j1[k] = MIN(jstep * j, ipixels_down - 1); j2[k] = MIN((jstep * (j+1)) - 1, ipixels_down - 1); k++; } } /* shuffle the work */ for (i=0; i<500; i++) { ii = rand() % 400; jj = rand() % 400; swap(&i1[ii], &i1[jj]); swap(&i2[ii], &i2[jj]); swap(&j1[ii], &j1[jj]); swap(&j2[ii], &j2[jj]); } /* send a piece of work to each worker (there must be more work than workers) */ k = 0; for (i=0; i<num_children; i++) { work[0] = k+1; work[1] = i1[k]; /* imin */ work[2] = i2[k]; /* imax */ work[3] = j1[k]; /* jmin */ work[4] = j2[k]; /* jmax */ /*printf("sending work(%d) to %d\n", k+1, i);fflush(stdout);*/ MPI_Send(work, 5, MPI_INT, 0, 100, child_comm[i]); MPI_Irecv(&result[i], 5, MPI_INT, 0, 200, child_comm[i], &child_request[i]); k++; } /* receive the results and hand out more work until the image is complete */ while (k<400) { MPI_Waitany(num_children, child_request, &index, &status); memcpy(work, &result[index], 5 * sizeof(int)); /*printf("master receiving data in k<400 loop.\n");fflush(stdout);*/ MPI_Recv(in_grid_array, (work[2] + 1 - work[1]) * (work[4] + 1 - work[3]), MPI_INT, 0, 201, child_comm[index], &status); /* draw data */ output_data(in_grid_array, &work[1], out_grid_array, ipixels_across, ipixels_down); work[0] = k+1; work[1] = i1[k]; work[2] = i2[k]; work[3] = j1[k]; work[4] = j2[k]; /*printf("sending work(%d) to %d\n", k+1, index);fflush(stdout);*/ MPI_Send(work, 5, MPI_INT, 0, 100, child_comm[index]); MPI_Irecv(&result[index], 5, MPI_INT, 0, 200, child_comm[index], &child_request[index]); k++; } /* receive the last pieces of work */ /* and tell everyone to stop */ for (i=0; i<num_children; i++) { MPI_Wait(&child_request[i], &status); memcpy(work, &result[i], 5 * sizeof(int)); /*printf("master receiving data in tail loop.\n");fflush(stdout);*/ MPI_Recv(in_grid_array, (work[2] + 1 - work[1]) * (work[4] + 1 - work[3]), MPI_INT, 0, 201, child_comm[i], &status); /* draw data */ output_data(in_grid_array, &work[1], out_grid_array, ipixels_across, ipixels_down); work[0] = 0; work[1] = 0; work[2] = 0; work[3] = 0; work[4] = 0; /*printf("sending %d to tell %d to stop\n", work[0], i);fflush(stdout);*/ MPI_Send(work, 5, MPI_INT, 0, 100, child_comm[i]); } /* tell the visualizer the image is done */ if (!use_stdin) { work[0] = 0; work[1] = 0; work[2] = 0; work[3] = 0; sock_write(sock, work, 4 * sizeof(int)); } } } else { for (;;) { /*printf("slave[%d] receiveing bounds.\n", pid);fflush(stdout);*/ MPI_Recv(&x_min, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&x_max, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&y_min, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&y_max, 1, MPI_DOUBLE, 0, 0, parent, MPI_STATUS_IGNORE); MPI_Recv(&imax_iterations, 1, MPI_INT, 0, 0, parent, MPI_STATUS_IGNORE); /*printf("slave[%d] received bounding box: (%f,%f)(%f,%f)\n", pid, x_min, y_min, x_max, y_max);fflush(stdout);*/ /* check for the end condition */ if (x_min == x_max && y_min == y_max) { /*printf("slave[%d] done.\n", pid);fflush(stdout);*/ break; } x_resolution = (x_max-x_min)/ ((double)ipixels_across); y_resolution = (y_max-y_min)/ ((double)ipixels_down); MPI_Recv(work, 5, MPI_INT, 0, 100, parent, &status); /*printf("slave[%d] received work: %d, (%d,%d)(%d,%d)\n", pid, work[0], work[1], work[2], work[3], work[4]);fflush(stdout);*/ while (work[0] != 0) { imin = work[1]; imax = work[2]; jmin = work[3]; jmax = work[4]; k = 0; for (j=jmin; j<=jmax; ++j) { coord_point.imaginary = y_max - j*y_resolution; /* go top to bottom */ for (i=imin; i<=imax; ++i) { /* Call Mandelbrot routine for each code, fill array with number of iterations. */ coord_point.real = x_min + i*x_resolution; /* go left to right */ if (julia == 1) { /* doing Julia set */ /* julia eq: z = z^2 + c, z_0 = grid coordinate, c = constant */ icount = single_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else if (alternate_equation == 1) { /* doing experimental form 1 */ icount = subtractive_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else if (alternate_equation == 2) { /* doing experimental form 2 */ icount = additive_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else { /* default to doing Mandelbrot set */ /* mandelbrot eq: z = z^2 + c, z_0 = c, c = grid coordinate */ icount = single_mandelbrot_point(coord_point, coord_point, imax_iterations, divergent_limit); } in_grid_array[k] = icount; ++k; } } /* send the result to the root */ /*printf("slave[%d] sending work %d back.\n", pid, work[0]);fflush(stdout);*/ MPI_Send(work, 5, MPI_INT, 0, 200, parent); /*printf("slave[%d] sending work %d data.\n", pid, work[0]);fflush(stdout);*/ MPI_Send(in_grid_array, (work[2] + 1 - work[1]) * (work[4] + 1 - work[3]), MPI_INT, 0, 201, parent); /* get the next piece of work */ /*printf("slave[%d] receiving new work.\n", pid);fflush(stdout);*/ MPI_Recv(work, 5, MPI_INT, 0, 100, parent, &status); /*printf("slave[%d] received work: %d, (%d,%d)(%d,%d)\n", pid, work[0], work[1], work[2], work[3], work[4]);fflush(stdout);*/ } } } if (master && save_image) { imax_iterations = 0; for (i=0; i<ipixels_across * ipixels_down; ++i) { /* look for "brightest" pixel value, for image use */ if (out_grid_array[i] > imax_iterations) imax_iterations = out_grid_array[i]; } if (julia == 0) printf("Done calculating mandelbrot, now creating file\n"); else printf("Done calculating julia, now creating file\n"); fflush(stdout); /* Print out the array in some appropriate form. */ if (julia == 0) { /* it's a mandelbrot */ sprintf(file_message, "Mandelbrot over (%lf-%lf,%lf-%lf), size %d x %d", x_min, x_max, y_min, y_max, ipixels_across, ipixels_down); } else { /* it's a julia */ sprintf(file_message, "Julia over (%lf-%lf,%lf-%lf), size %d x %d, center (%lf, %lf)", x_min, x_max, y_min, y_max, ipixels_across, ipixels_down, julia_constant.real, julia_constant.imaginary); } dumpimage(filename, out_grid_array, ipixels_across, ipixels_down, imax_iterations, file_message, num_colors, colors); } if (master) { for (i=0; i<num_children; i++) { MPI_Comm_disconnect(&child_comm[i]); } free(child_comm); free(child_request); free(colors); } MPI_Finalize(); return 0; }
int main(int argc, char *argv[]) { complex_t coord_point, julia_constant; double x_max, x_min, y_max, y_min, x_resolution, y_resolution; double divergent_limit; char file_message[160]; char filename[100]; int icount, imax_iterations; int ipixels_across, ipixels_down; int i, j, k, julia, alternate_equation; int imin, imax, jmin, jmax; int *work; /* make an integer array of size [N x M] to hold answers. */ int *grid_array = NULL; int numprocs; int namelen; char processor_name[MPI_MAX_PROCESSOR_NAME]; int num_colors; color_t *colors = NULL; int listener; int save_image = 0; int optval; int big_size; MPI_Win win; int error; int done; int use_datatypes = 1; MPI_Init(&argc, &argv); MPI_Comm_size(MPI_COMM_WORLD, &numprocs); MPI_Comm_rank(MPI_COMM_WORLD, &myid); MPI_Get_processor_name(processor_name, &namelen); if (numprocs == 1) { PrintUsage(); MPI_Finalize(); exit(0); } if (myid == 0) { printf("Welcome to the Mandelbrot/Julia set explorer.\n"); /* Get inputs-- region to view (must be within x/ymin to x/ymax, make sure xmax>xmin and ymax>ymin) and resolution (number of pixels along an edge, N x M, i.e. 256x256) */ read_mand_args(argc, argv, &imax_iterations, &ipixels_across, &ipixels_down, &x_min, &x_max, &y_min, &y_max, &julia, &julia_constant.real, &julia_constant.imaginary, &divergent_limit, &alternate_equation, filename, &num_colors, &use_stdin, &save_image, &use_datatypes); check_mand_params(&imax_iterations, &ipixels_across, &ipixels_down, &x_min, &x_max, &y_min, &y_max, &divergent_limit); if (julia == 1) /* we're doing a julia figure */ check_julia_params(&julia_constant.real, &julia_constant.imaginary); MPI_Bcast(&num_colors, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&imax_iterations, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&ipixels_across, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&ipixels_down, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&divergent_limit, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&julia, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&julia_constant.real, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&julia_constant.imaginary, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&alternate_equation, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&use_datatypes, 1, MPI_INT, 0, MPI_COMM_WORLD); } else { MPI_Bcast(&num_colors, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&imax_iterations, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&ipixels_across, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&ipixels_down, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&divergent_limit, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&julia, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&julia_constant.real, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&julia_constant.imaginary, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&alternate_equation, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&use_datatypes, 1, MPI_INT, 0, MPI_COMM_WORLD); } if (myid == 0) { colors = malloc((num_colors+1)* sizeof(color_t)); if (colors == NULL) { MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } Make_color_array(num_colors, colors); colors[num_colors] = 0; /* add one on the top to avoid edge errors */ } /* allocate memory */ big_size = ipixels_across * ipixels_down * sizeof(int); if (myid == 0) { /* window memory should be allocated by MPI in case the provider requires it */ error = MPI_Alloc_mem(big_size, MPI_INFO_NULL, &grid_array); if (error != MPI_SUCCESS) { printf("Memory allocation failed for data array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } /* allocate an array to put the workers tasks in */ work = (int*)malloc(numprocs * sizeof(int) * 5); if (work == NULL) { printf("Memory allocation failed for work array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } } else { /* the non-root processes just need scratch space to store data in */ if ( (grid_array = (int *)calloc(big_size, 1)) == NULL) { printf("Memory allocation failed for data array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } /* window memory should be allocated by MPI in case the provider requires it */ error = MPI_Alloc_mem(5 * sizeof(int), MPI_INFO_NULL, &work); if (error != MPI_SUCCESS) { printf("Memory allocation failed for work array, aborting.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } } if (myid == 0) { int istep, jstep; int i1[400], i2[400], j1[400], j2[400]; int ii, jj; struct sockaddr_in addr; int len; char line[1024], *token; srand(getpid()); if (!use_stdin) { addr.sin_family = AF_INET; addr.sin_addr.s_addr = INADDR_ANY; addr.sin_port = htons(DEFAULT_PORT); listener = socket(AF_INET, SOCK_STREAM, 0); if (listener == -1) { printf("unable to create a listener socket.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } if (bind(listener, &addr, sizeof(addr)) == -1) { addr.sin_port = 0; if (bind(listener, &addr, sizeof(addr)) == -1) { printf("unable to create a listener socket.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } } if (listen(listener, 1) == -1) { printf("unable to listen.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } len = sizeof(addr); getsockname(listener, &addr, &len); printf("%s listening on port %d\n", processor_name, ntohs(addr.sin_port)); fflush(stdout); sock = accept(listener, NULL, NULL); if (sock == -1) { printf("unable to accept a socket connection.\n"); MPI_Abort(MPI_COMM_WORLD, -1); exit(-1); } printf("accepted connection from visualization program.\n"); fflush(stdout); #ifdef HAVE_WINDOWS_H optval = 1; setsockopt(sock, IPPROTO_TCP, TCP_NODELAY, (char *)&optval, sizeof(optval)); #endif printf("sending image size to visualizer.\n"); sock_write(sock, &ipixels_across, sizeof(int)); sock_write(sock, &ipixels_down, sizeof(int)); sock_write(sock, &num_colors, sizeof(int)); sock_write(sock, &imax_iterations, sizeof(int)); } error = MPI_Win_create(grid_array, big_size, sizeof(int), MPI_INFO_NULL, MPI_COMM_WORLD, &win); if (error != MPI_SUCCESS) { printf("MPI_Win_create failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } for (;;) { /* get x_min, x_max, y_min, and y_max */ if (use_stdin) { printf("input xmin ymin xmax ymax max_iter, (0 0 0 0 0 to quit):\n");fflush(stdout); fgets(line, 1024, stdin); printf("read <%s> from stdin\n", line);fflush(stdout); token = strtok(line, " \n"); x_min = atof(token); token = strtok(NULL, " \n"); y_min = atof(token); token = strtok(NULL, " \n"); x_max = atof(token); token = strtok(NULL, " \n"); y_max = atof(token); token = strtok(NULL, " \n"); imax_iterations = atoi(token); } else { printf("reading xmin,ymin,xmax,ymax.\n");fflush(stdout); sock_read(sock, &x_min, sizeof(double)); sock_read(sock, &y_min, sizeof(double)); sock_read(sock, &x_max, sizeof(double)); sock_read(sock, &y_max, sizeof(double)); sock_read(sock, &imax_iterations, sizeof(int)); } printf("x0,y0 = (%f, %f) x1,y1 = (%f,%f) max_iter = %d\n", x_min, y_min, x_max, y_max, imax_iterations);fflush(stdout); /* break the work up into 400 pieces */ istep = ipixels_across / 20; jstep = ipixels_down / 20; if (istep < 1) istep = 1; if (jstep < 1) jstep = 1; k = 0; for (i=0; i<20; i++) { for (j=0; j<20; j++) { i1[k] = MIN(istep * i, ipixels_across - 1); i2[k] = MIN((istep * (i+1)) - 1, ipixels_across - 1); j1[k] = MIN(jstep * j, ipixels_down - 1); j2[k] = MIN((jstep * (j+1)) - 1, ipixels_down - 1); k++; } } /* shuffle the work */ for (i=0; i<500; i++) { ii = rand() % 400; jj = rand() % 400; swap(&i1[ii], &i1[jj]); swap(&i2[ii], &i2[jj]); swap(&j1[ii], &j1[jj]); swap(&j2[ii], &j2[jj]); } /*printf("bcasting the limits: (%f,%f)(%f,%f)\n", x_min, y_min, x_max, y_max);fflush(stdout);*/ /* let everyone know the limits */ MPI_Bcast(&x_min, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&x_max, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&y_min, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&y_max, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&imax_iterations, 1, MPI_INT, 0, MPI_COMM_WORLD); /* check for the end condition */ if (x_min == x_max && y_min == y_max) { break; } /* put one piece of work to each worker for each epoch until the work is exhausted */ k = 0; done = 0; while (!done) { error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'handout work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* hand out work */ for (i=1; i<numprocs; i++) { if (!done) { work[(i*5)+0] = k+1; work[(i*5)+1] = i1[k]; /* imin */ work[(i*5)+2] = i2[k]; /* imax */ work[(i*5)+3] = j1[k]; /* jmin */ work[(i*5)+4] = j2[k]; /* jmax */ } else { work[(i*5)+0] = -1; work[(i*5)+1] = -1; work[(i*5)+2] = -1; work[(i*5)+3] = -1; work[(i*5)+4] = -1; } /*printf("sending work(%d) to %d\n", k+1, cur_proc);fflush(stdout);*/ error = MPI_Put(&work[i*5], 5, MPI_INT, i, 0, 5, MPI_INT, win); if (error != MPI_SUCCESS) { printf("put failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } if (k<399) k++; else done = 1; } error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'handout work' -> 'do work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* do work */ error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'do work' -> 'collect results' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* send the results to the visualizer */ for (i=1; i<numprocs; i++) { if (work[i*5] != -1) { sock_write(sock, &work[i*5 + 1], 4*sizeof(int)); for (j=work[i*5+3]; j<=work[i*5+4]; j++) { sock_write(sock, &grid_array[(j*ipixels_across)+work[i*5+1]], (work[i*5+2]+1-work[i*5+1])*sizeof(int)); } } } } error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'collect results' -> 'done work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* hand out "done" work */ for (i=1; i<numprocs; i++) { work[(i*5)+0] = 0; work[(i*5)+1] = 0; work[(i*5)+2] = 0; work[(i*5)+3] = 0; work[(i*5)+4] = 0; error = MPI_Put(&work[i*5], 5, MPI_INT, i, 0, 5, MPI_INT, win); if (error != MPI_SUCCESS) { printf("put failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } } error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'done work' -> 'done' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* tell the visualizer the image is done */ if (!use_stdin) { work[0] = 0; work[1] = 0; work[2] = 0; work[3] = 0; sock_write(sock, work, 4 * sizeof(int)); } } } else { MPI_Datatype dtype; error = MPI_Win_create(work, 5*sizeof(int), sizeof(int), MPI_INFO_NULL, MPI_COMM_WORLD, &win); if (error != MPI_SUCCESS) { printf("MPI_Win_create failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } for (;;) { MPI_Bcast(&x_min, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&x_max, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&y_min, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&y_max, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&imax_iterations, 1, MPI_INT, 0, MPI_COMM_WORLD); /* check for the end condition */ if (x_min == x_max && y_min == y_max) { break; } x_resolution = (x_max-x_min)/ ((double)ipixels_across); y_resolution = (y_max-y_min)/ ((double)ipixels_down); error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'receive work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* receive work from the root */ error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'receive work' -> 'do work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } while (work[0] != 0) { imin = work[1]; imax = work[2]; jmin = work[3]; jmax = work[4]; if (use_datatypes) { MPI_Type_vector(jmax - jmin + 1, /* rows */ imax - imin + 1, /* column width */ ipixels_across, /* stride, distance between rows */ MPI_INT, &dtype); MPI_Type_commit(&dtype); k = 0; } for (j=jmin; j<=jmax; ++j) { coord_point.imaginary = y_max - j*y_resolution; /* go top to bottom */ for (i=imin; i<=imax; ++i) { /* Call Mandelbrot routine for each code, fill array with number of iterations. */ coord_point.real = x_min + i*x_resolution; /* go left to right */ if (julia == 1) { /* doing Julia set */ /* julia eq: z = z^2 + c, z_0 = grid coordinate, c = constant */ icount = single_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else if (alternate_equation == 1) { /* doing experimental form 1 */ icount = subtractive_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else if (alternate_equation == 2) { /* doing experimental form 2 */ icount = additive_mandelbrot_point(coord_point, julia_constant, imax_iterations, divergent_limit); } else { /* default to doing Mandelbrot set */ /* mandelbrot eq: z = z^2 + c, z_0 = c, c = grid coordinate */ icount = single_mandelbrot_point(coord_point, coord_point, imax_iterations, divergent_limit); } if (use_datatypes) { grid_array[k++] = icount; } else { grid_array[(j*ipixels_across) + i] = icount; error = MPI_Put(&grid_array[(j*ipixels_across) + i], 1, MPI_INT, 0, (j * ipixels_across) + i, 1, MPI_INT, win); if (error != MPI_SUCCESS) { printf("put failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } } } } if (use_datatypes) { MPI_Put(grid_array, k, MPI_INT, 0, (jmin * ipixels_across) + imin, 1, dtype, win); } /* synch with the root */ error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'do work' -> 'wait for work to be collected' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } if (use_datatypes) { MPI_Type_free(&dtype); } /* fence while the root writes to the visualizer. */ error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'wait for work to be collected' -> 'receive work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } /* fence to allow the root to put the next piece of work */ error = MPI_Win_fence(0, win); if (error != MPI_SUCCESS) { printf("'receive work' -> 'do work' fence failed, error 0x%x\n", error); MPI_Abort(MPI_COMM_WORLD, -1); } } } } if (myid == 0 && save_image) { imax_iterations = 0; for (i=0; i<ipixels_across * ipixels_down; ++i) { /* look for "brightest" pixel value, for image use */ if (grid_array[i] > imax_iterations) imax_iterations = grid_array[i]; } if (julia == 0) printf("Done calculating mandelbrot, now creating file\n"); else printf("Done calculating julia, now creating file\n"); fflush(stdout); /* Print out the array in some appropriate form. */ if (julia == 0) { /* it's a mandelbrot */ sprintf(file_message, "Mandelbrot over (%lf-%lf,%lf-%lf), size %d x %d", x_min, x_max, y_min, y_max, ipixels_across, ipixels_down); } else { /* it's a julia */ sprintf(file_message, "Julia over (%lf-%lf,%lf-%lf), size %d x %d, center (%lf, %lf)", x_min, x_max, y_min, y_max, ipixels_across, ipixels_down, julia_constant.real, julia_constant.imaginary); } dumpimage(filename, grid_array, ipixels_across, ipixels_down, imax_iterations, file_message, num_colors, colors); } MPI_Finalize(); if (colors) free(colors); return 0; } /* end of main */