File: qmr.sci

package info (click to toggle)
scilab 5.2.2-9
  • links: PTS, VCS
  • area: main
  • in suites: squeeze
  • size: 334,832 kB
  • ctags: 52,586
  • sloc: xml: 526,945; ansic: 223,590; fortran: 163,080; java: 56,934; cpp: 33,840; tcl: 27,936; sh: 20,397; makefile: 9,908; ml: 9,451; perl: 1,323; cs: 614; lisp: 30
file content (393 lines) | stat: -rw-r--r-- 10,697 bytes parent folder | download
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
// Scilab ( http://www.scilab.org/ ) - This file is part of Scilab
// Copyright (C) XXXX-2008 - INRIA
// Copyright (C) 2005 - IRISA - Sage Group

// 
// This file must be used under the terms of the CeCILL.
// This source file is licensed as described in the file COPYING, which
// you should have received as part of this distribution.  The terms
// are also available at    
// http://www.cecill.info/licences/Licence_CeCILL_V2-en.txt


// [x, flag, err, iter, res] = qmr( A, Ap, b, x, M1, M2, max_it, tol )
//
// QMR solves the linear system Ax=b using the 
// Quasi Minimal Residual method with preconditioning.
//
// input   A        REAL matrix or function
//         x        REAL initial guess vector
//         b        REAL right hand side vector
//         M1       REAL left preconditioner matrix
//         M2       REAL right preconditioner matrix
//         max_it   INTEGER maximum number of iterations
//         tol      REAL error tolerance
//
// output  x        REAL solution vector
//         flag     INTEGER: 0: solution found to tolerance
//                           1: no convergence given max_it
//                     breakdown:
//                          -1: rho
//                          -2: Beta
//                          -3: gam
//                          -4: delta
//                          -5: ep
//                          -6: xi
//         err      REAL final residual norm
//         iter     INTEGER number of iterations performed
//         res      REAL residual vector

//     Details of this algorithm are described in 
//
//     "Templates for the Solution of Linear Systems: Building Blocks 
//     for Iterative Methods", 
//     Barrett, Berry, Chan, Demmel, Donato, Dongarra, Eijkhout,
//     Pozo, Romine, and Van der Vorst, SIAM Publications, 1993
//     (ftp netlib2.cs.utk.edu; cd linalg; get templates.ps).
//
//     "Iterative Methods for Sparse Linear Systems, Second Edition"
//     Saad, SIAM Publications, 2003
//     (ftp ftp.cs.umn.edu; cd dept/users/saad/PS; get all_ps.zip).

// Sage Group (IRISA, 2005)

function [x, flag, err, iter, res] = qmr( A, varargin)

// -----------------------
// Parsing input arguments
// -----------------------
  [lhs,rhs]=argn(0);
  if ( rhs < 2 ),
    error(msprintf(gettext("%s: Wrong number of input arguments: At least %d expected.\n"),"qmr",2));
  end

  // Parsing the matrix A
  select type(A)
  case 1 then
    cpt=1;
  case 5 then
    cpt=1;
  case 13 then
    cpt=0;
  else
    error(msprintf(gettext("%s: Wrong type for input argument #%d.\n"),"qmr",1));
  end

  // If A is a matrix (dense or sparse)
  if (cpt==1),
    if (size(A,1) ~= size(A,2)),
      error(msprintf(gettext("%s: Wrong type for input argument #%d: Square matrix expected.\n"),"qmr",1));
    end
    deff('y=matvec(x)','y=A*x');
    deff('y=matvecp(x)','y=A''*x');
    fct=0;
  end

  // If A is a function
  if (cpt==0),
    if (rhs == 1),
      error(msprintf(gettext("%s: Wrong type for input argument #%d: Transpose of the function %s expected.\n"),"qmr",1,"A"));
    end
    matvec=A;
    fct=1;
  end
  if (rhs >= 2 & fct==1 ),
    matvecp=varargin(1);
    if ( type(matvecp) ~= 13 ),
      error(msprintf(gettext("%s: Wrong value for input argument #%d: Transpose of the function %s expected.\n"),"qmr",2,"A"));
    end
  end

  // Parsing right hand side b
  if ( rhs >= fct+2 ),
    b=varargin(fct+1);
    if ( size(b,2) ~= 1),
      error(msprintf(gettext("%s: Wrong value for input argument #%d: Column vector expected.\n"),"qmr",3));
    end
  else 
    error(msprintf(gettext("%s: Wrong value for input argument #%d: Vector expected.\n"),"qmr",3));
  end

  // Parsing initial vector x
  if ( rhs >= fct+3),
    x=varargin(fct+2);
    if (size(x,2) ~= 1),
      error(msprintf(gettext("%s: Wrong type for input argument #%d: Column vector expected.\n"),"qmr",4));
    end
    if ( size(x,1) ~= size(b,1)),
      error(msprintf(gettext("%s: Wrong type for input argument #%d: Same size as the input argument #%d expected.\n"),"qmr",4,3));
    end
  else
    x=zeros(size(b,1),1);
  end

  //--------------------------------------------------------
  // Parsing of the preconditioner matrix M1
  //--------------------------------------------------------

  if (rhs >=fct+4),
    Prec_g=varargin(fct+3);
    select type(Prec_g)
    case 1 then
      cpt=1;
    case 5 then
      cpt=1;
    case 13 then
      cpt=0;
    end 
    if ( cpt==1 ),
	  if (size(Prec_g,1) ~= size(Prec_g,2)),
		error(msprintf(gettext("%s: Wrong type for input argument #%d: Square matrix expected.\n"),"qmr",5));
      end 
      if (size(Prec_g,1)~=size(b,1)),
		error(msprintf(gettext("%s: Wrong size for input argument #%d: Must be same size as input argument #%d"),"qmr",5,3));
      end 
      deff('y=precond_g(x)','y=Prec_g \ x');
      deff('y=precondp_g(x)','y=Prec_g'' \ x');
    end
    if ( cpt==0 ),
      if ( rhs >= fct+5 ),
	Precp_g=varargin(fct+4);
	select type(Precp_g)
	case 1 then
	  cpt1=1;
	case 5 then
	  cpt1=1;
	case 13 then
	  cpt1=0;
	end 
	if ( cpt1==0 ),
	  precond_g=Prec_g;
	  precondp_g=Precp_g;
	  fct=fct+1;
	else
	  error(msprintf(gettext("%s: Wrong type for input argument #%d: Transpose of the function %s expected.\n"),"qmr",5,"M1"));
	end
  else
	error(msprintf(gettext("%s: Wrong type for input argument #%d: Transpose of the function %s expected.\n"),"qmr",5,"M1"));
  end
    end
  else
    deff('y=precond_g(x)','y=x');
    deff('y=precondp_g(x)','y=x');
  end  
  
  //--------------------------------------------------------
  // Parsing of the preconditioner matrix M1
  //--------------------------------------------------------

  if (rhs >=fct+5),
    Prec_d=varargin(fct+4);
    select type(Prec_d)
    case 1 then
      cpt=1;
    case 5 then
      cpt=1;
    case 13 then
      cpt=0;
    end 
    if ( cpt==1 ),
      if (size(Prec_d,1) ~= size(Prec_d,2)),
		error(msprintf(gettext("%s: Wrong type for input argument #%d: Square matrix expected.\n"),"qmr",6));
      end 
      if (size(Prec_d,1)~=size(b,1)),
		error(msprintf(gettext("%s: Wrong size for input argument #%d: Must be same size as input argument #%d"),"qmr",6,3));
      end 
      deff('y=precond_d(x)','y=Prec_d \ x');
      deff('y=precondp_d(x)','y=Prec_d'' \ x');
    end
    if ( cpt==0 ),
      if ( rhs >= fct+6 ),
	Precp_d=varargin(fct+5);
	select type(Precp_d)
	case 1 then
	  cpt1=1;
	case 5 then
	  cpt1=1;
	case 13 then
	  cpt1=0;
	end 
	if ( cpt1==0 ),
	  precond_d=Prec_d;
	  precondp_d=Precp_d;
	  fct=fct+1;
	else
	  error(msprintf(gettext("%s: Wrong type for input argument #%d: Transpose of the function %s expected.\n"),"qmr",6,"M2"));

	end
  else
		  error(msprintf(gettext("%s: Wrong type for input argument #%d: Transpose of the function %s expected.\n"),"qmr",6,"M2"));
      end
    end
    
  else
    deff('y=precond_d(x)','y=x');
    deff('y=precondp_d(x)','y=x');
  end

  //--------------------------------------------------------
  // Parsing of the maximum number of iterations max_it
  //--------------------------------------------------------

  if (rhs >= fct+6),
    max_it=varargin(fct+5);
    if (size(max_it,1) ~= 1 | size(max_it,2) ~=1),
	  error(msprintf(gettext("%s: Wrong type for input argument #%d: Scalar expected.\n"),"qmr",7));

    end 
  else
    max_it=size(b,1);
  end

  //--------------------------------------------------------
  // Parsing of the error tolerance tol
  //--------------------------------------------------------

  if (rhs == fct+7),
    tol=varargin(fct+6);
    if (size(tol,1) ~= 1 | size(tol,2) ~=1),
	  error(msprintf(gettext("%s: Wrong type for input argument #%d: Scalar expected.\n"),"qmr",8));
	  
    end
  else
    tol=1000*%eps;
  end

  //--------------------------------------------------------
  // test about input arguments number
  //--------------------------------------------------------

  if (rhs > fct+8),
	error(msprintf(gettext("%s: Wrong number of input arguments: %d to %d expected.\n"),"qmr",2,8));
  end

  // ------------
  // Computations
  // ------------

  // initialization
  i = 0;
  flag = 0;
  bnrm2 = norm( b );
  if  (bnrm2 == 0.0), 
    bnrm2 = 1.0; 
  end

  //   r = b - A*x;
  r = b - matvec(x);
  err = norm( r ) / bnrm2;
  res = err;
  if ( err < tol ), return; end

  // [M1,M2] = lu( M );

  v_tld = r;
  // y = M1 \ v_tld;
  y = precond_g(v_tld)
  rho = norm( y );

  w_tld = r;
  //   z = M2' \ w_tld;
  z = precondp_d(w_tld);
  xi = norm( z );

  gam =  1.0;
  eta = -1.0;
  theta =  0.0;

  for i = 1:max_it,                      // begin iteration

    if ( rho == 0.0 | xi == 0.0 ), iter=i; break; end

    v = v_tld / rho;
    y = y / rho;

    w = w_tld / xi;
    z = z / xi;

    delta = z'*y;
    if ( delta == 0.0 ), iter=i; break; end

    //    y_tld = M2 \ y;
    y_tld = precond_d(y);
    //    z_tld = M1'\ z;
    z_tld = precondp_g(z);

    if ( i > 1 ),                       // direction vector 
      p = y_tld - ( xi*delta / ep )*p;
      q = z_tld - ( rho*delta / ep )*q;
    else
      p = y_tld;
      q = z_tld;
    end

    //    p_tld = A*p;
    p_tld = matvec(p);

    ep = q'*p_tld;
    if ( ep == 0.0 ), iter=i; break; end

    Beta = ep / delta;
    if ( Beta == 0.0 ), iter=i; break; end

    v_tld = p_tld - Beta*v;
    //    y =  M1 \ v_tld;
    y = precond_g(v_tld);
    
    rho_1 = rho;
    rho = norm( y );
    //    w_tld = ( A'*q ) - ( Beta*w );
    w_tld = ( matvecp(q) ) - ( Beta*w );
    //    z =  M2' \ w_tld;
    z =  precondp_d(w_tld);

    xi = norm( z );

    gamma_1 = gam;
    theta_1 = theta;

    theta = rho / ( gamma_1*Beta );
    gam = 1.0 / sqrt( 1.0 + (theta^2) );
    if ( gam == 0.0 ), iter=i; break; end

    eta = -eta*rho_1*(gam^2) / ( Beta*(gamma_1^2) );

    if ( i > 1 ),                         // compute adjustment
      d = eta*p + (( theta_1*gam )^2)*d;
      s = eta*p_tld + (( theta_1*gam )^2)*s;
    else
      d = eta*p;
      s = eta*p_tld;
    end

    x = x + d;                               // update approximation

    r = r - s;                               // update residual
    err = norm( r ) / bnrm2;               // check convergence
    res = [res;err];
    
    if ( err <= tol ), iter=i; break; end

    if ( i == max_it ), iter=i; end
    
  end

  if ( err <= tol ),                        // converged
    flag =  0;
  elseif ( rho == 0.0 ),                      // breakdown
    flag = -1;
  elseif ( Beta == 0.0 ),
    flag = -2;
  elseif ( gam == 0.0 ),
    flag = -3;
  elseif ( delta == 0.0 ),
    flag = -4;
  elseif ( ep == 0.0 ),
    flag = -5;
  elseif ( xi == 0.0 ),
    flag = -6;
  else                                        // no convergence
    flag = 1;
  end

endfunction