File: spgmr.c

package info (click to toggle)
xppaut 8.0%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 7,332 kB
  • sloc: ansic: 74,690; makefile: 127; sh: 92
file content (443 lines) | stat: -rw-r--r-- 11,568 bytes parent folder | download | duplicates (3)
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
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
/******************************************************************
 *                                                                *
 * File          : spgmr.c                                        *
 * Programmers   : Scott D. Cohen and Alan C. Hindmarsh @ LLNL    *
 * Last Modified : 1 September 1994                               *
 *----------------------------------------------------------------*
 * This is the implementation file for the scaled preconditioned  *
 * GMRES (SPGMR) iterative linear solver.                         *
 *                                                                *
 ******************************************************************/


#include <stdio.h>
#include <stdlib.h>
#include "iterativ.h"
#include "spgmr.h"
#include "llnltyps.h"
#include "vector.h"
#include "llnlmath.h"


#define ZERO RCONST(0.0)
#define ONE  RCONST(1.0)


/*************** Private Helper Function Prototype *******************/

static void FreeVectorArray(N_Vector *A, int indMax);
 

/* Implementation of Spgmr algorithm */


/*************** SpgmrMalloc *****************************************/

SpgmrMem SpgmrMalloc(integer N, int l_max, void *machEnv)
{
  SpgmrMem mem;
  N_Vector *V, xcor, vtemp;
  real **Hes, *givens, *yg;
  int k, i;
 
  /* Check the input parameters */

  if ((N <= 0) || (l_max <= 0)) return(NULL);

  /* Get memory for the Krylov basis vectors V[0], ..., V[l_max] */
  
  V = (N_Vector *) malloc((l_max+1)*sizeof(N_Vector));
  if (V == NULL) return(NULL);

  for (k=0; k <= l_max; k++) {
    V[k] = N_VNew(N, machEnv);
    if (V[k] == NULL) {
      FreeVectorArray(V, k-1);
      return(NULL);
    }
  }

  /* Get memory for the Hessenberg matrix Hes */

  Hes = (real **) malloc((l_max+1)*sizeof(real *)); 
  if (Hes == NULL) {
    FreeVectorArray(V, l_max);
    return(NULL);
  }

  for (k=0; k <= l_max; k++) {
    Hes[k] = (real *) malloc(l_max*sizeof(real));
    if (Hes[k] == NULL) {
      for (i=0; i < k; i++) free(Hes[i]);
      FreeVectorArray(V, l_max);
      return(NULL);
    }
  }

  /* Get memory for Givens rotation components */

  givens = (real *) malloc(2*l_max*sizeof(real));
  if (givens == NULL) {
    for (i=0; i <= l_max; i++) free(Hes[i]);
    FreeVectorArray(V, l_max);
    return(NULL);
  }

  /* Get memory to hold the correction to z_tilde */

  xcor = N_VNew(N, machEnv);
  if (xcor == NULL) {
    free(givens);
    for (i=0; i <= l_max; i++) free(Hes[i]);
    FreeVectorArray(V, l_max);
    return(NULL);
  }

  /* Get memory to hold SPGMR y and g vectors */

  yg = (real *) malloc((l_max+1)*sizeof(real));
  if (yg == NULL) {
    N_VFree(xcor);
    free(givens);
    for (i=0; i <= l_max; i++) free(Hes[i]);
    FreeVectorArray(V, l_max);
    return(NULL);
  }

  /* Get an array to hold a temporary vector */

  vtemp = N_VNew(N, machEnv);
  if (vtemp == NULL) {
    free(yg);
    N_VFree(xcor);
    free(givens);
    for (i=0; i <= l_max; i++) free(Hes[i]);
    FreeVectorArray(V, l_max);
    return(NULL);
  }

  /* Get memory for an SpgmrMemRec containing SPGMR matrices and vectors */

  mem = (SpgmrMem) malloc(sizeof(SpgmrMemRec));
  if (mem == NULL) {
    N_VFree(vtemp);
    free(yg);
    N_VFree(xcor);
    free(givens);
    for (i=0; i <= l_max; i++) free(Hes[i]);
    FreeVectorArray(V, l_max);
    return(NULL); 
  }

  /* Set the fields of mem */

  mem->N = N;
  mem->l_max = l_max;
  mem->V = V;
  mem->Hes = Hes;
  mem->givens = givens;
  mem->xcor = xcor;
  mem->yg = yg;
  mem->vtemp = vtemp;

  /* Return the pointer to SPGMR memory */

  return(mem);
}


/*************** SpgmrSolve ******************************************/

int SpgmrSolve(SpgmrMem mem, void *A_data, N_Vector x, N_Vector b,
               int pretype, int gstype, real delta, int max_restarts,
	       void *P_data, N_Vector sx, N_Vector sb, ATimesFn atimes,
	       PSolveFn psolve, real *res_norm, int *nli, int *nps)
{
  N_Vector *V, xcor, vtemp;
  real **Hes, *givens, *yg;
  /*real s_r0_norm, beta, rotation_product, r_norm, s_product, rho;*/
  real  beta, rotation_product, r_norm, s_product, rho=0.0;
  bool preOnLeft, preOnRight, scale_x, scale_b, converged;
  int i, j, k, l, l_plus_1, l_max, krydim=0, ier, ntries;
  
  if (mem == NULL) return(SPGMR_MEM_NULL);
  
  /* Make local copies of mem variables */

  l_max  = mem->l_max;
  V      = mem->V;
  Hes    = mem->Hes;
  givens = mem->givens;
  xcor   = mem->xcor;
  yg     = mem->yg;
  vtemp  = mem->vtemp;
  
  *nli = *nps = 0;     /* Initialize counters */
  converged = FALSE;   /* Initialize converged flag */
  
  if (max_restarts < 0) max_restarts = 0;
  
  if ((pretype != LEFT) && (pretype != RIGHT) && (pretype != BOTH))
    pretype = NONE;
  
  preOnLeft  = ((pretype == LEFT) || (pretype == BOTH));
  preOnRight = ((pretype == RIGHT) || (pretype == BOTH));
  scale_x    = (sx != NULL);
  scale_b    = (sb != NULL);

  /* Set vtemp and V[0] to initial (unscaled) residual r_0 = b - A*x_0  */

  if (N_VDotProd(x, x) == ZERO) {
    N_VScale(ONE, b, vtemp);  
  } else {
    if (atimes(A_data, x, vtemp) != 0)   
      return(SPGMR_ATIMES_FAIL);
    N_VLinearSum(ONE, b, -ONE, vtemp, vtemp);  
  }
  N_VScale(ONE, vtemp, V[0]);

  /* Apply b-scaling to vtemp, get L2 norm of sb r_0, and return if small */
/*
  if (scale_b) N_VProd(sb, vtemp, vtemp);
  s_r0_norm = RSqrt(N_VDotProd(vtemp, vtemp));
  if (s_r0_norm <= delta) return(SPGMR_SUCCESS);
*/  
  /* Apply left preconditioner and b-scaling to V[0] = r_0 */
  
  if (preOnLeft) {
    ier = psolve(P_data, V[0], vtemp, LEFT);
    (*nps)++;
    if (ier != 0)
      return((ier < 0) ? SPGMR_PSOLVE_FAIL_UNREC : SPGMR_PSOLVE_FAIL_REC);
  } else {
    N_VScale(ONE, V[0], vtemp);
  }
  
  if (scale_b) {
    N_VProd(sb, vtemp, V[0]);   
  } else {
    N_VScale(ONE, vtemp, V[0]);
  }

  /* Set r_norm = beta to L2 norm of V[0] = sb P1_inv r_0, and
     return if small  */

  *res_norm = r_norm = beta = RSqrt(N_VDotProd(V[0], V[0])); 
  if (r_norm <= delta)
    return(SPGMR_SUCCESS);

  /* Set xcor = 0 */

  N_VConst(ZERO, xcor);


  /* Begin outer iterations: up to (max_restarts + 1) attempts */
  
  for (ntries = 0; ntries <= max_restarts; ntries++) {

    /* Initialize the Hessenberg matrix Hes and Givens rotation
       product.  Normalize the initial vector V[0].             */
   
    for (i=0; i <= l_max; i++)
      for (j=0; j < l_max; j++)
	Hes[i][j] = ZERO;

    rotation_product = ONE;
    
    N_VScale(ONE/r_norm, V[0], V[0]);

    /* Inner loop: generate Krylov sequence and Arnoldi basis */
    
    for(l=0; l < l_max; l++) {

      (*nli)++;

      krydim = l_plus_1 = l + 1;
      
      /* Generate A-tilde V[l], where A-tilde = sb P1_inv A P2_inv sx_inv */

      /* Apply x-scaling: vtemp = sx_inv V[l] */
      if (scale_x) {
	N_VDiv(V[l], sx, vtemp);
      } else {
	N_VScale(ONE, V[l], vtemp);
      }

      /* Apply right precoditioner: vtemp = P2_inv sx_inv V[l] */ 
      N_VScale(ONE, vtemp, V[l_plus_1]);
      if (preOnRight) {
	ier = psolve(P_data, V[l_plus_1], vtemp, RIGHT);
	(*nps)++;
	if (ier != 0)
	  return((ier < 0) ? SPGMR_PSOLVE_FAIL_UNREC : SPGMR_PSOLVE_FAIL_REC);
      }

      /* Apply A: V[l+1] = A P2_inv sx_inv V[l] */      
      if (atimes(A_data, vtemp, V[l_plus_1] ) != 0) 
	return(SPGMR_ATIMES_FAIL);

      /* Apply left preconditioning: vtemp = P1_inv A P2_inv sx_inv V[l] */
      if (preOnLeft) {
	ier = psolve(P_data, V[l_plus_1], vtemp, LEFT);
	(*nps)++;
	if (ier != 0)
	  return((ier < 0) ? SPGMR_PSOLVE_FAIL_UNREC : SPGMR_PSOLVE_FAIL_REC);
      } else {
	N_VScale(ONE, V[l_plus_1], vtemp);
      }

      /* Apply b-scaling: V[l+1] = sb P1_inv A P2_inv sx_inv V[l] */
      if (scale_b) {
	N_VProd(sb, vtemp, V[l_plus_1]);
      } else {
	N_VScale(ONE, vtemp, V[l_plus_1]);
      }
      
      /*  Orthogonalize V[l+1] against previous V[i]: V[l+1] = w_tilde. */

      if (gstype == CLASSICAL_GS) {
	if (ClassicalGS(V, Hes, l_plus_1, l_max, &(Hes[l_plus_1][l]),
			vtemp, yg) != 0)
	  return(SPGMR_GS_FAIL);
      } else {
	if (ModifiedGS(V, Hes, l_plus_1, l_max, &(Hes[l_plus_1][l])) != 0) 
	  return(SPGMR_GS_FAIL);
      }

      /*  Update the QR factorization of Hes  */

      if(QRfact(krydim, Hes, givens, l) != 0 )
	return(SPGMR_QRFACT_FAIL);

      /*  Update residual norm estimate; break if convergence test passes */
      
      rotation_product *= givens[2*l+1];
    
      if ((*res_norm = rho = ABS(rotation_product*r_norm)) <= delta) {
	converged = TRUE;
	break;
      }
      
      /* Normalize V[l+1] with norm value from the Gram-Schmidt */
      N_VScale(ONE/Hes[l_plus_1][l], V[l_plus_1], V[l_plus_1]);
    }
    
    /* Inner loop is done.  Compute the new correction vector xcor */

    /* Construct g, then solve for y */
    yg[0] = r_norm;
    for (i=1; i <= krydim; i++) yg[i]=ZERO;
    if (QRsol(krydim, Hes, givens, yg) != 0)
      return(SPGMR_QRSOL_FAIL);
    
    /* Add correction vector V_l y to xcor */
    for (k=0; k < krydim; k++)
      N_VLinearSum(yg[k], V[k], ONE, xcor, xcor);

    /* If converged, construct the final solution vector x */
    if (converged) {
     
      /* Apply x-scaling and right precond.: vtemp = P2_inv sx_inv xcor */
  
      if (scale_x) N_VDiv(xcor, sx, xcor);
      if (preOnRight) {
	ier = psolve(P_data, xcor, vtemp, RIGHT);
	(*nps)++;
	if (ier != 0)
	   return((ier < 0) ? SPGMR_PSOLVE_FAIL_UNREC : SPGMR_PSOLVE_FAIL_REC);
      } else {
	N_VScale(ONE, xcor, vtemp);
      }

      /* Add correction to initial x to get final solution x, and return */
      N_VLinearSum(ONE, x, ONE, vtemp, x);

      return(SPGMR_SUCCESS);
    }

    /* Not yet converged; if allowed, prepare for restart */

    if (ntries == max_restarts) break;

    /* Construct last column of Q in yg */
    s_product = ONE;
    for (i=krydim; i > 0; i--) {
      yg[i] = s_product*givens[2*i-2];
      s_product *= givens[2*i-1];
    }
    yg[0] = s_product;

    /* Scale r_norm and yg */
    r_norm *= s_product;
    for (i=0; i <= krydim; i++)
      yg[i] *= r_norm;
    r_norm = ABS(r_norm);

    /* Multiply yg by V_(krydim+1) to get last residual vector; restart */
    N_VScale(yg[0], V[0], V[0]);
    for( k=1; k <= krydim; k++)
      N_VLinearSum(yg[k], V[k], ONE, V[0], V[0]);

  }

  /* Failed to converge, even after allowed restarts.
     If the residual norm was reduced below its initial value, compute
     and return x anyway.  Otherwise return failure flag.              */

  if (rho < beta) {

    /* Apply the x-scaling and right precond.: vtemp = P2_inv sx_inv xcor */
    
    if (scale_x) N_VDiv(xcor, sx, xcor);
    if (preOnRight) {
      ier = psolve(P_data, xcor, vtemp, RIGHT);
      (*nps)++;
      if (ier != 0)
	return((ier < 0) ? SPGMR_PSOLVE_FAIL_UNREC : SPGMR_PSOLVE_FAIL_REC);
    } else {
      N_VScale(ONE, xcor, vtemp);
    }

    /* Add vtemp to initial x to get final solution x, and return */
    N_VLinearSum(ONE, x, ONE, vtemp, x);
    
    return(SPGMR_RES_REDUCED);
  }

  return(SPGMR_CONV_FAIL); 
}

/*************** SpgmrFree *******************************************/

void SpgmrFree(SpgmrMem mem)
{
  int i, l_max;
  real **Hes;

  if (mem == NULL) return;

  l_max = mem->l_max;
  Hes = mem->Hes;

  FreeVectorArray(mem->V, l_max);
  for (i=0; i <= l_max; i++) free(Hes[i]);
  free(Hes);
  free(mem->givens);
  N_VFree(mem->xcor);
  free(mem->yg);
  N_VFree(mem->vtemp);

  free(mem);
}


/*************** Private Helper Function: FreeVectorArray ************/

static void FreeVectorArray(N_Vector *A, int indMax)
{
  int j;

  for (j=0; j <= indMax; j++) N_VFree(A[j]);

  free(A);
}