File: SparseCmplxLU.cc

package info (click to toggle)
octave3.0 1%3A3.0.1-6lenny3
  • links: PTS, VCS
  • area: main
  • in suites: lenny
  • size: 49,912 kB
  • ctags: 27,149
  • sloc: cpp: 166,567; fortran: 61,399; ansic: 8,766; sh: 5,856; lex: 4,419; makefile: 3,965; yacc: 3,013; lisp: 1,692; objc: 889; perl: 795; awk: 532
file content (445 lines) | stat: -rw-r--r-- 12,130 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
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
444
445
/*

Copyright (C) 2004, 2005, 2006, 2007 David Bateman
Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004 Andy Adler

This file is part of Octave.

Octave is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 3 of the License, or (at your
option) any later version.

Octave is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
for more details.

You should have received a copy of the GNU General Public License
along with Octave; see the file COPYING.  If not, see
<http://www.gnu.org/licenses/>.

*/

#ifdef HAVE_CONFIG_H
#include <config.h>
#endif

#include <vector>

#include "lo-error.h"

#include "SparseCmplxLU.h"
#include "oct-spparms.h"

// Instantiate the base LU class for the types we need.

#include "sparse-base-lu.h"
#include "sparse-base-lu.cc"

template class sparse_base_lu <SparseComplexMatrix, Complex, SparseMatrix, double>;

#include "oct-sparse.h"

SparseComplexLU::SparseComplexLU (const SparseComplexMatrix& a, 
				  double piv_thres)
{
#ifdef HAVE_UMFPACK
  octave_idx_type nr = a.rows ();
  octave_idx_type nc = a.cols ();

  // Setup the control parameters
  Matrix Control (UMFPACK_CONTROL, 1);
  double *control = Control.fortran_vec ();
  UMFPACK_ZNAME (defaults) (control);

  double tmp = octave_sparse_params::get_key ("spumoni");
  if (!xisnan (tmp))
    Control (UMFPACK_PRL) = tmp;
  if (piv_thres >= 0.)
    {
      piv_thres = (piv_thres > 1. ? 1. : piv_thres);
      Control (UMFPACK_SYM_PIVOT_TOLERANCE) = piv_thres;
      Control (UMFPACK_PIVOT_TOLERANCE) = piv_thres;
    }
  else
    {
      tmp = octave_sparse_params::get_key ("piv_tol");
      if (!xisnan (tmp))
	{
	  Control (UMFPACK_SYM_PIVOT_TOLERANCE) = tmp;
      Control (UMFPACK_PIVOT_TOLERANCE) = tmp;
	}
    }

  // Set whether we are allowed to modify Q or not
  tmp = octave_sparse_params::get_key ("autoamd");
  if (!xisnan (tmp))
    Control (UMFPACK_FIXQ) = tmp;

  // Turn-off UMFPACK scaling for LU 
  Control (UMFPACK_SCALE) = UMFPACK_SCALE_NONE;

  UMFPACK_ZNAME (report_control) (control);

  const octave_idx_type *Ap = a.cidx ();
  const octave_idx_type *Ai = a.ridx ();
  const Complex *Ax = a.data ();

  UMFPACK_ZNAME (report_matrix) (nr, nc, Ap, Ai,
				 reinterpret_cast<const double *> (Ax),
				 NULL, 1, control);

  void *Symbolic;
  Matrix Info (1, UMFPACK_INFO);
  double *info = Info.fortran_vec ();
  int status = UMFPACK_ZNAME (qsymbolic) (nr, nc, Ap, Ai, 
					  reinterpret_cast<const double *> (Ax),
					  NULL, NULL,
					  &Symbolic, control, info);

  if (status < 0)
    {
      (*current_liboctave_error_handler) 
	    ("SparseComplexLU::SparseComplexLU symbolic factorization failed");

      UMFPACK_ZNAME (report_status) (control, status);
      UMFPACK_ZNAME (report_info) (control, info);

      UMFPACK_ZNAME (free_symbolic) (&Symbolic) ;
    }
  else
    {
      UMFPACK_ZNAME (report_symbolic) (Symbolic, control);

      void *Numeric;
      status = UMFPACK_ZNAME (numeric) (Ap, Ai,
					reinterpret_cast<const double *> (Ax),
					NULL, Symbolic, &Numeric, control,
					info);
      UMFPACK_ZNAME (free_symbolic) (&Symbolic) ;

      cond = Info (UMFPACK_RCOND);

      if (status < 0)
	{
	  (*current_liboctave_error_handler) 
	    ("SparseComplexLU::SparseComplexLU numeric factorization failed");

	  UMFPACK_ZNAME (report_status) (control, status);
	  UMFPACK_ZNAME (report_info) (control, info);

	  UMFPACK_ZNAME (free_numeric) (&Numeric);
	}
      else
	{
	  UMFPACK_ZNAME (report_numeric) (Numeric, control);

	  octave_idx_type lnz, unz, ignore1, ignore2, ignore3;
	  status = UMFPACK_ZNAME (get_lunz) (&lnz, &unz, &ignore1,
					&ignore2, &ignore3, Numeric) ;
	  
	  if (status < 0)
	    {
	      (*current_liboctave_error_handler) 
		("SparseComplexLU::SparseComplexLU extracting LU factors failed");

	      UMFPACK_ZNAME (report_status) (control, status);
	      UMFPACK_ZNAME (report_info) (control, info);

	      UMFPACK_ZNAME (free_numeric) (&Numeric);
	    }
	  else
	    {
	      octave_idx_type n_inner = (nr < nc ? nr : nc);

	      if (lnz < 1)
		Lfact = SparseComplexMatrix (n_inner, nr,
					     static_cast<octave_idx_type> (1));
	      else
		Lfact = SparseComplexMatrix (n_inner, nr, lnz);

	      octave_idx_type *Ltp = Lfact.cidx ();
	      octave_idx_type *Ltj = Lfact.ridx ();
	      Complex *Ltx = Lfact.data ();

	      if (unz < 1)
		Ufact = SparseComplexMatrix (n_inner, nc,
					     static_cast<octave_idx_type> (1));
	      else
		Ufact = SparseComplexMatrix (n_inner, nc, unz);

	      octave_idx_type *Up = Ufact.cidx ();
	      octave_idx_type *Uj = Ufact.ridx ();
	      Complex *Ux = Ufact.data ();
	      
	      P.resize (nr);
	      octave_idx_type *p = P.fortran_vec ();

	      Q.resize (nc);
	      octave_idx_type *q = Q.fortran_vec ();

	      octave_idx_type do_recip;
	      status = UMFPACK_ZNAME (get_numeric) (Ltp, Ltj,
						    reinterpret_cast<double *> (Ltx),
						    NULL, Up, Uj,
						    reinterpret_cast <double *> (Ux),
						    NULL, p, q, NULL, NULL,
						    &do_recip, NULL, Numeric);

	      UMFPACK_ZNAME (free_numeric) (&Numeric) ;

	      if (status < 0 || do_recip)
		{
		  (*current_liboctave_error_handler) 
		    ("SparseComplexLU::SparseComplexLU extracting LU factors failed");

		  UMFPACK_ZNAME (report_status) (control, status);
		}
	      else
		{
		  Lfact = Lfact.transpose ();

		  UMFPACK_ZNAME (report_matrix) (nr, n_inner,
					    Lfact.cidx (), Lfact.ridx (), 
					    reinterpret_cast<double *> (Lfact.data()), 
					    NULL, 1, control);

		  UMFPACK_ZNAME (report_matrix) (n_inner, nc,
					    Ufact.cidx (), Ufact.ridx (), 
					    reinterpret_cast<double *> (Ufact.data()), 
					    NULL, 1, control);
		  UMFPACK_ZNAME (report_perm) (nr, p, control);
		  UMFPACK_ZNAME (report_perm) (nc, q, control);
		}

	      UMFPACK_ZNAME (report_info) (control, info);
	    }
	}
    }
#else
  (*current_liboctave_error_handler) ("UMFPACK not installed");
#endif
}

SparseComplexLU::SparseComplexLU (const SparseComplexMatrix& a, 
				  const ColumnVector& Qinit, 
				  double piv_thres, bool FixedQ,
				  double droptol, bool milu, bool udiag)
{
#ifdef HAVE_UMFPACK
  if (milu)
    (*current_liboctave_error_handler) 
      ("Modified incomplete LU not implemented");   
  else
    {
      octave_idx_type nr = a.rows ();
      octave_idx_type nc = a.cols ();

      // Setup the control parameters
      Matrix Control (UMFPACK_CONTROL, 1);
      double *control = Control.fortran_vec ();
      UMFPACK_ZNAME (defaults) (control);

      double tmp = octave_sparse_params::get_key ("spumoni");
      if (!xisnan (tmp))
	Control (UMFPACK_PRL) = tmp;
      if (piv_thres >= 0.)
	{
	  piv_thres = (piv_thres > 1. ? 1. : piv_thres);
	  Control (UMFPACK_SYM_PIVOT_TOLERANCE) = piv_thres;
	  Control (UMFPACK_PIVOT_TOLERANCE) = piv_thres;
	}
      else
	{
	  tmp = octave_sparse_params::get_key ("piv_tol");
	  if (!xisnan (tmp))
	    {
	      Control (UMFPACK_SYM_PIVOT_TOLERANCE) = tmp;
	      Control (UMFPACK_PIVOT_TOLERANCE) = tmp;
	    }
	}

      if (droptol >= 0.)
	Control (UMFPACK_DROPTOL) = droptol;

      // Set whether we are allowed to modify Q or not
      if (FixedQ)
	Control (UMFPACK_FIXQ) = 1.0;
      else
	{
	  tmp = octave_sparse_params::get_key ("autoamd");
	  if (!xisnan (tmp))
	    Control (UMFPACK_FIXQ) = tmp;
	}

      // Turn-off UMFPACK scaling for LU 
      Control (UMFPACK_SCALE) = UMFPACK_SCALE_NONE;

      UMFPACK_ZNAME (report_control) (control);

      const octave_idx_type *Ap = a.cidx ();
      const octave_idx_type *Ai = a.ridx ();
      const Complex *Ax = a.data ();

      UMFPACK_ZNAME (report_matrix) (nr, nc, Ap, Ai, 
				reinterpret_cast<const double *> (Ax), NULL,
				1, control);

      void *Symbolic;
      Matrix Info (1, UMFPACK_INFO);
      double *info = Info.fortran_vec ();
      int status;

      // Null loop so that qinit is imediately deallocated when not
      // needed
      do {
	OCTAVE_LOCAL_BUFFER (octave_idx_type, qinit, nc);

	for (octave_idx_type i = 0; i < nc; i++)
	  qinit [i] = static_cast<octave_idx_type> (Qinit (i));

	status = UMFPACK_ZNAME (qsymbolic) (nr, nc, Ap, Ai, 
				       reinterpret_cast<const double *> (Ax),
				       NULL, qinit, &Symbolic, control, 
				       info);
      } while (0);

      if (status < 0)
	{
	  (*current_liboctave_error_handler) 
	    ("SparseComplexLU::SparseComplexLU symbolic factorization failed");

	  UMFPACK_ZNAME (report_status) (control, status);
	  UMFPACK_ZNAME (report_info) (control, info);

	  UMFPACK_ZNAME (free_symbolic) (&Symbolic) ;
	}
      else
	{
	  UMFPACK_ZNAME (report_symbolic) (Symbolic, control);

	  void *Numeric;
	  status = UMFPACK_ZNAME (numeric) (Ap, Ai, 
				       reinterpret_cast<const double *> (Ax), NULL,
				       Symbolic, &Numeric, control, info) ;
	  UMFPACK_ZNAME (free_symbolic) (&Symbolic) ;

	  cond = Info (UMFPACK_RCOND);

	  if (status < 0)
	    {
	      (*current_liboctave_error_handler) 
		("SparseComplexLU::SparseComplexLU numeric factorization failed");

	      UMFPACK_ZNAME (report_status) (control, status);
	      UMFPACK_ZNAME (report_info) (control, info);

	      UMFPACK_ZNAME (free_numeric) (&Numeric);
	    }
	  else
	    {
	      UMFPACK_ZNAME (report_numeric) (Numeric, control);

	      octave_idx_type lnz, unz, ignore1, ignore2, ignore3;
	      status = UMFPACK_ZNAME (get_lunz) (&lnz, &unz,
					    &ignore1, &ignore2, &ignore3, Numeric);
	  
	      if (status < 0)
		{
		  (*current_liboctave_error_handler) 
		    ("SparseComplexLU::SparseComplexLU extracting LU factors failed");

		  UMFPACK_ZNAME (report_status) (control, status);
		  UMFPACK_ZNAME (report_info) (control, info);

		  UMFPACK_ZNAME (free_numeric) (&Numeric);
		}
	      else
		{
		  octave_idx_type n_inner = (nr < nc ? nr : nc);

		  if (lnz < 1)
		    Lfact = SparseComplexMatrix (n_inner, nr,
		       static_cast<octave_idx_type> (1));
		  else
		    Lfact = SparseComplexMatrix (n_inner, nr, lnz);

		  octave_idx_type *Ltp = Lfact.cidx ();
		  octave_idx_type *Ltj = Lfact.ridx ();
		  Complex *Ltx = Lfact.data ();

		  if (unz < 1)
		    Ufact = SparseComplexMatrix (n_inner, nc,
		       static_cast<octave_idx_type> (1));
		  else
		    Ufact = SparseComplexMatrix  (n_inner, nc, unz);

		  octave_idx_type *Up = Ufact.cidx ();
		  octave_idx_type *Uj = Ufact.ridx ();
		  Complex *Ux = Ufact.data ();
	      
		  P.resize (nr);
		  octave_idx_type *p = P.fortran_vec ();

		  Q.resize (nc);
		  octave_idx_type *q = Q.fortran_vec ();

		  octave_idx_type do_recip;
		  status = 
		    UMFPACK_ZNAME (get_numeric) (Ltp, Ltj, 
					    reinterpret_cast<double *> (Ltx),
					    NULL, Up, Uj,
					    reinterpret_cast<double *> (Ux), 
					    NULL, p, q, NULL, NULL, 
					    &do_recip, NULL, Numeric) ;

		  UMFPACK_ZNAME (free_numeric) (&Numeric) ;

		  if (status < 0 || do_recip)
		    {
		      (*current_liboctave_error_handler) 
			("SparseComplexLU::SparseComplexLU extracting LU factors failed");

		      UMFPACK_ZNAME (report_status) (control, 
								     status);
		    }
		  else
		    {
		      Lfact = Lfact.transpose ();

		      UMFPACK_ZNAME (report_matrix) (nr, n_inner, 
						Lfact.cidx (), 
						Lfact.ridx (), 
						reinterpret_cast<double *> (Lfact.data()), 
						NULL, 1, control);

		      UMFPACK_ZNAME (report_matrix) (n_inner, nc, 
						Ufact.cidx (), 
						Ufact.ridx (), 
						reinterpret_cast<double *> (Ufact.data()), 
						NULL, 1, control);
		      UMFPACK_ZNAME (report_perm) (nr, p, control);
		      UMFPACK_ZNAME (report_perm) (nc, q, control);
		    }

		  UMFPACK_ZNAME (report_info) (control, info);
		}
	    }
	}

      if (udiag)
	(*current_liboctave_error_handler) 
	  ("Option udiag of incomplete LU not implemented");   
    }
#else
  (*current_liboctave_error_handler) ("UMFPACK not installed");
#endif
}

/*
;;; Local Variables: ***
;;; mode: C++ ***
;;; End: ***
*/