1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
/* -------------------------------------------------------
 * Solve a distributed lasso problem, i.e.,
 *
 *   minimize \lambda * ||x||_1 + 0.5 * ||Ax - b||_2^2
 *
 * The implementation uses MPI for distributed communication
 * and the GNU Scientific Library (GSL) for math.
 * Compile: make
 * run code: mpiexec -n 1 ./IST data_directory
 *
 * Author:   Zhimin Peng
 * Date:     01/11/2013
 * Modified: 02/06/2013
 *--------------------------------------------------------*/

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include "mmio.h"
#include <mpi.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_linalg.h>

void shrink(gsl_vector *v, double sigma);
double objective(gsl_vector *x, double lambda, gsl_vector *z);

int main(int argc, char **argv) {
  
  const int MAX_ITER      = 50; // number of iteration
  const double TOL        = 1e-4; // tolerence
  const double lambda_tol = 1e-4; // tolerence for update lambda
  int rank; // process ID
  int size; // number of processes

  MPI_Init(&argc, &argv); // initialize MPI environment
  MPI_Comm_rank(MPI_COMM_WORLD, &rank); // determine current running process
  MPI_Comm_size(MPI_COMM_WORLD, &size); // total number of processes


  char* dataCenterDir = "Your data directory";
  char* big_dir; // directory of data
  if(argc==2)
    big_dir = argv[1];
  else
    big_dir = "big1";

  /* Read in local data */
  FILE *f, *test;
  int m, n, row, col;
  double entry, startTime, endTime, commStartTime, commEndTime, commTime;
  char s[100];

  /* -------------------------------------------------------
   * Subsystem n will look for files called An.dat and bn.dat
   * in the current directory; these are its local data and do 
   * not need to be visible to any other processes. Note that
   * m and n here refer to the dimensions of the 
   * local coefficient matrix.
   * -------------------------------------------------------*/

  /* Read A */
  sprintf(s, "%s/%s/A%d.dat",dataCenterDir,big_dir, rank + 1);
  printf("[%d] reading %s\n", rank, s);  
  f = fopen(s, "r");
  if (f == NULL) {
    printf("[%d] ERROR: %s does not exist, exiting.\n", rank, s);
   exit(EXIT_FAILURE);
  }
  mm_read_mtx_array_size(f, &m, &n);
  gsl_matrix *A = gsl_matrix_calloc(m, n);
  for (int i = 0; i < m*n; i++) {
    row = i % m;
    col = floor(i/m);
    fscanf(f, "%lf", &entry);
    gsl_matrix_set(A, row, col, entry);
  }
  fclose(f);

  /* Read b */
  sprintf(s, "%s/%s/b.dat", dataCenterDir, big_dir);
  f = fopen(s, "r");
  if (f == NULL) {
    printf("[%d] ERROR: %s does not exist, exiting.\n", rank, s);
    exit(EXIT_FAILURE);
  }
  mm_read_mtx_array_size(f, &m, &n);
  gsl_vector *b = gsl_vector_calloc(m);
  for (int i = 0; i < m; i++) {
    fscanf(f, "%lf", &entry);
    gsl_vector_set(b, i, entry);
  }
  fclose(f);
  
  /* Read xs */
  sprintf(s, "%s/%s/xs%d.dat", dataCenterDir, big_dir, rank + 1);
  printf("[%d] reading %s\n", rank, s);
  f = fopen(s, "r");
  if (f == NULL) {
    printf("[%d] ERROR: %s does not exist, exiting.\n", rank, s);
    exit(EXIT_FAILURE);
  }
  mm_read_mtx_array_size(f, &m, &n);
  gsl_vector *xs = gsl_vector_calloc(m);
  for (int i = 0; i < m; i++) {
    fscanf(f, "%lf", &entry);
    gsl_vector_set(xs, i, entry);
  }
  fclose(f);
  // [m, n] = size(A);
  m = A->size1;
  n = A->size2;
  MPI_Barrier(MPI_COMM_WORLD);
 
  // These are all variables related to FISTA
  gsl_vector *x      = gsl_vector_calloc(n);
  gsl_vector *u      = gsl_vector_calloc(n);
  gsl_vector *xold   = gsl_vector_calloc(n);
  gsl_vector *z      = gsl_vector_calloc(m);
  gsl_vector *y      = gsl_vector_calloc(m);
  gsl_vector *w      = gsl_vector_calloc(n);
  gsl_vector *Ax     = gsl_vector_calloc(m);
  gsl_vector *xdiff  = gsl_vector_calloc(n);
  double send[1]; // an array used to aggregate 3 scalars at once
  double recv[1]; // used to receive the results of these aggregations
  double err;
  double xs_local_nrm[1], xs_nrm[1]; // calculate the two norm of xs
  xs_local_nrm[0] = gsl_blas_dnrm2(xs);
  xs_local_nrm[0] = xs_local_nrm[0]* xs_local_nrm[0];
  MPI_Allreduce(xs_local_nrm, xs_nrm, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
  xs_nrm[0] = sqrt(xs_nrm[0]);

  // FISTA parameters
  double delta, t1=1, t2=1;
  double lambda = 0.01;

  // approximate ||A||_2 by power method
  double sum_x = 0, err0=0, err1;
  double tol= 1e-6,local_normx, normx, normy, local_err1;
  int cnt=0;
  // calcuate x = sum(abs(A), 1)';
  startTime = MPI_Wtime();
   for(int j=0;j<n;j++){
    gsl_matrix_get_col(z, A, j);
    sum_x = gsl_blas_dasum(z);
    gsl_vector_set(x, j, sum_x);
  }
  local_err1 = gsl_blas_dnrm2(x);
  local_err1 = local_err1*local_err1;
  MPI_Allreduce(&local_err1, &err1,  1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);  // y = sum A * x
  err1 = sqrt(err1);
  gsl_vector_scale(x, 1.0/err1);
  while(fabs(err1 - err0) > tol * err1){
    err0 = err1;
    gsl_blas_dgemv(CblasNoTrans, 1, A, x, 0, Ax); // Ax = A*x
    MPI_Allreduce(Ax->data, y->data,  m, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); // y = sum A * x
    gsl_blas_dgemv(CblasTrans, 1, A, y, 0, x); 
    local_normx = gsl_blas_dnrm2(x);
    local_normx = local_normx * local_normx;
    MPI_Allreduce(&local_normx, &normx,  1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); // y = sum A * x
    normx = sqrt(normx);
    normy = gsl_blas_dnrm2(y);
  
    err1 = normx/normy;
    gsl_vector_scale(x, 1.0/normx);
    cnt = cnt+1;
    if(cnt > 100){
      break;
    }
  }
  endTime = MPI_Wtime();
  printf("spertral norm evaluation time: %e \n", endTime - startTime);
  
  /*---------------------- 
     initialize variables
    ----------------------*/
  gsl_vector_set_zero(x);
  gsl_vector_set_zero(y);
  gsl_vector_set_zero(z);
  
  delta = 1.00/(err1*err1);

  if (rank == 0) {
    printf("%3s %10s %10s %10s %10s %10s\n", "#", "r norm", "eps_pri", "s norm", "eps_dual", "objective");
    printf("delta = %lf \n", delta);

    sprintf(s, "Results/test%d.m", size);
    test = fopen(s, "w");
    fprintf(test,"res = [ \n");
  }
  
  int iter = 0;
  if (rank == 0) {
    printf("%3s %10s %10s\n", "#", "err", "objective");
  }

  /* Main FISTA solver loop */
  while (iter < MAX_ITER) {

    startTime = MPI_Wtime();

    t1 = t2;
    gsl_vector_memcpy(xold, x); // copy x to old x;
    gsl_blas_dgemv(CblasNoTrans, 1, A, u, 0, Ax); // A_i x_i = A_i * x_i

    commStartTime = MPI_Wtime();
    MPI_Allreduce(Ax->data, z->data,  m, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); // z = A*x
    commEndTime = MPI_Wtime();
    commTime = commEndTime - commStartTime; // calculate the communication time
    gsl_vector_sub(z, b);
    gsl_blas_dgemv(CblasTrans, 1, A, z, 0, w); // w = A' (Ax - b)   
    gsl_vector_scale(w, delta);  // w = delta * w
    gsl_vector_sub(u, w);  // u = x - delta * A'(Ax - b)
    gsl_vector_memcpy(x, u); 
    shrink(x, delta * lambda); // shrink(x, alpha*lambda)
 
    // FISTA 
    t2 = 0.5 + 0.5*sqrt(1+4*t1*t1);
    gsl_vector_sub(xold, x);
    err = gsl_blas_dnrm2(xold);
    err = err*err;
    gsl_vector_scale(xold, (1 - t1)/t2);
    gsl_vector_memcpy(u, x);
    gsl_vector_add(u, xold);

    /* check the error */
    gsl_vector_memcpy(xdiff,xs);
    gsl_vector_sub(xdiff, x);
    err = gsl_blas_dnrm2(xdiff);
    send[0] = err*err;
    MPI_Allreduce(send, recv, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); 
    recv[0] = sqrt(recv[0])/xs_nrm[0];
  

    /* termination check */
    if (rank == 0) {
      endTime = MPI_Wtime();
      printf("%3d %e %10.4f %e \n", iter,
	     recv[0],  objective(x, lambda, z), lambda);
      fprintf(test, "%e %e %e %e;\n", recv[0], objective(x, lambda, z), 
	      endTime - startTime, commTime);
    }
    
    if(recv[0] < TOL){
      break;
    }

    iter++;
  }
  
  /* Have the master write out the results to disk */
  if (rank == 0) {
    fprintf(test,"] \n");
    fprintf(test,"totalTime = sum(res(:,3)) \n");
    fprintf(test,"commTime = sum(res(:,4)) \n");
    fclose(test);
    f = fopen("Results/solution.dat", "w");
    fprintf(f,"x = [ \n");
    gsl_vector_fprintf(f, x, "%lf");
    fprintf(f,"] \n");
    fclose(f);
 
    printf("Elapsed time is: %lf \n", endTime - startTime);
  }
  
  MPI_Finalize(); /* Shut down the MPI execution environment */
  
  /* Clear memory */
  gsl_matrix_free(A);
  gsl_vector_free(b);
  gsl_vector_free(x);
  gsl_vector_free(y);
  gsl_vector_free(z);
  gsl_vector_free(w);
  gsl_vector_free(xdiff);
  gsl_vector_free(Ax);
  return 0;
}

/* calculate objective function */
double objective(gsl_vector *x, double lambda, gsl_vector *z) {
  double obj = 0;
  double temp =0.0;
  temp = gsl_blas_dnrm2(z);
  temp = temp*temp/2;
  double foo;
  foo = gsl_blas_dasum(x);
  obj = lambda*foo + temp;
  return obj;
}

/* shrinkage function */
void shrink(gsl_vector *v, double sigma) {
  double vi;
  for (int i = 0; i < v->size; i++) {
    vi = gsl_vector_get(v, i);
    if (vi > sigma)       { gsl_vector_set(v, i, vi-sigma); }
    else if (vi < -sigma) { gsl_vector_set(v, i, vi+sigma); }
    else              { gsl_vector_set(v, i, 0); }
  }
}