MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
solvers.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.googlecode.com.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "linalg.hpp"
13 #include <iostream>
14 #include <iomanip>
15 #include <algorithm>
16 #include <cmath>
17 
18 namespace mfem
19 {
20 
21 using namespace std;
22 
24  : Solver(0, true)
25 {
26  oper = NULL;
27  prec = NULL;
28  max_iter = 10;
29  print_level = -1;
30  rel_tol = abs_tol = 0.0;
31 #ifdef MFEM_USE_MPI
32  dot_prod_type = 0;
33 #endif
34 }
35 
36 #ifdef MFEM_USE_MPI
38  : Solver(0, true)
39 {
40  oper = NULL;
41  prec = NULL;
42  max_iter = 10;
43  print_level = -1;
44  rel_tol = abs_tol = 0.0;
45  dot_prod_type = 1;
46  comm = _comm;
47 }
48 #endif
49 
50 double IterativeSolver::Dot(const Vector &x, const Vector &y) const
51 {
52 #ifndef MFEM_USE_MPI
53  return (x * y);
54 #else
55  if (dot_prod_type == 0)
56  {
57  return (x * y);
58  }
59  else
60  {
61  double local_dot = (x * y);
62  double global_dot;
63 
64  MPI_Allreduce(&local_dot, &global_dot, 1, MPI_DOUBLE, MPI_SUM, comm);
65 
66  return global_dot;
67  }
68 #endif
69 }
70 
71 void IterativeSolver::SetPrintLevel(int print_lvl)
72 {
73 #ifndef MFEM_USE_MPI
74  print_level = print_lvl;
75 #else
76  if (dot_prod_type == 0)
77  {
78  print_level = print_lvl;
79  }
80  else
81  {
82  int rank;
83  MPI_Comm_rank(comm, &rank);
84  if (rank == 0)
85  print_level = print_lvl;
86  }
87 #endif
88 }
89 
91 {
92  prec = &pr;
93  prec->iterative_mode = false;
94 }
95 
97 {
98  oper = &op;
99  height = op.Height();
100  width = op.Width();
101  if (prec)
102  prec->SetOperator(*oper);
103 }
104 
105 
107 {
108  r.SetSize(width);
109  z.SetSize(width);
110 }
111 
112 void SLISolver::Mult(const Vector &b, Vector &x) const
113 {
114  int i;
115 
116  // Optimized preconditioned SLI with fixed number of iterations and given
117  // initial guess
118  if (!rel_tol && iterative_mode && prec)
119  {
120  for (i = 0; i < max_iter; i++)
121  {
122  oper->Mult(x, r); // r = A x
123  subtract(b, r, r); // r = b - A x
124  prec->Mult(r, z); // z = B r
125  add(x, 1.0, z, x); // x = x + B (b - A x)
126  }
127  converged = 1;
128  final_iter = i;
129  return;
130  }
131 
132  // Optimized preconditioned SLI with fixed number of iterations and zero
133  // initial guess
134  if (!rel_tol && !iterative_mode && prec)
135  {
136  prec->Mult(b, x); // x = B b (initial guess 0)
137  for (i = 1; i < max_iter; i++)
138  {
139  oper->Mult(x, r); // r = A x
140  subtract(b, r, r); // r = b - A x
141  prec->Mult(r, z); // z = B r
142  add(x, 1.0, z, x); // x = x + B (b - A x)
143  }
144  converged = 1;
145  final_iter = i;
146  return;
147  }
148 
149  // General version of SLI with a relative tolerance, optional preconditioner
150  // and optional initial guess
151  double r0, nom, nom0, nomold = 1, cf;
152 
153  if (iterative_mode)
154  {
155  oper->Mult(x, r);
156  subtract(b, r, r); // r = b - A x
157  }
158  else
159  {
160  r = b;
161  x = 0.0;
162  }
163 
164  if (prec)
165  {
166  prec->Mult(r, z); // z = B r
167  nom0 = nom = Dot(z, r);
168  }
169  else
170  {
171  nom0 = nom = Dot(r, r);
172  }
173 
174  if (print_level == 1)
175  cout << " Iteration : " << setw(3) << 0 << " (B r, r) = "
176  << nom << '\n';
177 
178  r0 = fmax(nom*rel_tol*rel_tol, abs_tol*abs_tol);
179  if (nom <= r0)
180  {
181  converged = 1;
182  final_iter = 0;
183  final_norm = sqrt(nom);
184  return;
185  }
186 
187  // start iteration
188  converged = 0;
190  for (i = 1; true; )
191  {
192  if (prec) // x = x + B (b - A x)
193  add(x, 1.0, z, x);
194  else
195  add(x, 1.0, r, x);
196 
197  oper->Mult(x, r);
198  subtract(b, r, r); // r = b - A x
199 
200  if (prec)
201  {
202  prec->Mult(r, z); // z = B r
203  nom = Dot(z, r);
204  }
205  else
206  {
207  nom = Dot(r, r);
208  }
209 
210  cf = sqrt(nom/nomold);
211  if (print_level == 1)
212  cout << " Iteration : " << setw(3) << i << " (B r, r) = "
213  << nom << "\tConv. rate: " << cf << '\n';
214  nomold = nom;
215 
216  if (nom < r0)
217  {
218  if (print_level == 2)
219  cout << "Number of SLI iterations: " << i << '\n'
220  << "Conv. rate: " << cf << '\n';
221  else
222  if (print_level == 3)
223  cout << "(B r_0, r_0) = " << nom0 << '\n'
224  << "(B r_N, r_N) = " << nom << '\n'
225  << "Number of SLI iterations: " << i << '\n';
226  converged = 1;
227  final_iter = i;
228  break;
229  }
230 
231  if (++i > max_iter)
232  break;
233  }
234 
235  if (print_level >= 0 && !converged)
236  {
237  cerr << "SLI: No convergence!" << '\n';
238  cout << "(B r_0, r_0) = " << nom0 << '\n'
239  << "(B r_N, r_N) = " << nom << '\n'
240  << "Number of SLI iterations: " << final_iter << '\n';
241  }
242  if (print_level >= 1 || (print_level >= 0 && !converged))
243  {
244  cout << "Average reduction factor = "
245  << pow (nom/nom0, 0.5/final_iter) << '\n';
246  }
247  final_norm = sqrt(nom);
248 }
249 
250 void SLI(const Operator &A, const Vector &b, Vector &x,
251  int print_iter, int max_num_iter,
252  double RTOLERANCE, double ATOLERANCE)
253 {
254  SLISolver sli;
255  sli.SetPrintLevel(print_iter);
256  sli.SetMaxIter(max_num_iter);
257  sli.SetRelTol(sqrt(RTOLERANCE));
258  sli.SetAbsTol(sqrt(ATOLERANCE));
259  sli.SetOperator(A);
260  sli.Mult(b, x);
261 }
262 
263 void SLI(const Operator &A, Solver &B, const Vector &b, Vector &x,
264  int print_iter, int max_num_iter,
265  double RTOLERANCE, double ATOLERANCE)
266 {
267  SLISolver sli;
268  sli.SetPrintLevel(print_iter);
269  sli.SetMaxIter(max_num_iter);
270  sli.SetRelTol(sqrt(RTOLERANCE));
271  sli.SetAbsTol(sqrt(ATOLERANCE));
272  sli.SetOperator(A);
273  sli.SetPreconditioner(B);
274  sli.Mult(b, x);
275 }
276 
277 
279 {
280  r.SetSize(width);
281  d.SetSize(width);
282  z.SetSize(width);
283 }
284 
285 void CGSolver::Mult(const Vector &b, Vector &x) const
286 {
287  int i;
288  double r0, den, nom, nom0, betanom, alpha, beta;
289 
290  if (iterative_mode)
291  {
292  oper->Mult(x, r);
293  subtract(b, r, r); // r = b - A x
294  }
295  else
296  {
297  r = b;
298  x = 0.0;
299  }
300 
301  if (prec)
302  {
303  prec->Mult(r, z); // z = B r
304  d = z;
305  }
306  else
307  {
308  d = r;
309  }
310  nom0 = nom = Dot(d, r);
311 
312  if (print_level == 1)
313  cout << " Iteration : " << setw(3) << 0 << " (B r, r) = "
314  << nom << '\n';
315 
316  r0 = fmax(nom*rel_tol*rel_tol, abs_tol*abs_tol);
317  if (nom <= r0)
318  {
319  converged = 1;
320  final_iter = 0;
321  final_norm = sqrt(nom);
322  return;
323  }
324 
325  oper->Mult(d, z); // z = A d
326  den = Dot(z, d);
327 
328  if (print_level >= 0 && den < 0.0)
329  {
330  cout << "Negative denominator in step 0 of PCG: " << den << '\n';
331  }
332 
333  if (den == 0.0)
334  {
335  converged = 0;
336  final_iter = 0;
337  final_norm = sqrt(nom);
338  return;
339  }
340 
341  // start iteration
342  converged = 0;
344  for (i = 1; true; )
345  {
346  alpha = nom/den;
347  add(x, alpha, d, x); // x = x + alpha d
348  add(r, -alpha, z, r); // r = r - alpha A d
349 
350  if (prec)
351  {
352  prec->Mult(r, z); // z = B r
353  betanom = Dot(r, z);
354  }
355  else
356  {
357  betanom = Dot(r, r);
358  }
359 
360  if (print_level == 1)
361  cout << " Iteration : " << setw(3) << i << " (B r, r) = "
362  << betanom << '\n';
363 
364  if (betanom < r0)
365  {
366  if (print_level == 2)
367  cout << "Number of PCG iterations: " << i << '\n';
368  else
369  if (print_level == 3)
370  cout << "(B r_0, r_0) = " << nom0 << '\n'
371  << "(B r_N, r_N) = " << betanom << '\n'
372  << "Number of PCG iterations: " << i << '\n';
373  converged = 1;
374  final_iter = i;
375  break;
376  }
377 
378  if (++i > max_iter)
379  break;
380 
381  beta = betanom/nom;
382  if (prec)
383  add(z, beta, d, d); // d = z + beta d
384  else
385  add(r, beta, d, d);
386  oper->Mult(d, z); // z = A d
387  den = Dot(d, z);
388  if (den <= 0.0)
389  {
390  if (print_level >= 0 && Dot(d, d) > 0.0)
391  cout <<"PCG: The operator is not postive definite. (Ad, d) = "
392  << den << '\n';
393  }
394  nom = betanom;
395  }
396  if (print_level >= 0 && !converged)
397  {
398  cerr << "PCG: No convergence!" << '\n';
399  cout << "(B r_0, r_0) = " << nom0 << '\n'
400  << "(B r_N, r_N) = " << betanom << '\n'
401  << "Number of PCG iterations: " << final_iter << '\n';
402  }
403  if (print_level >= 1 || (print_level >= 0 && !converged))
404  {
405  cout << "Average reduction factor = "
406  << pow (betanom/nom0, 0.5/final_iter) << '\n';
407  }
408  final_norm = sqrt(betanom);
409 }
410 
411 void CG(const Operator &A, const Vector &b, Vector &x,
412  int print_iter, int max_num_iter,
413  double RTOLERANCE, double ATOLERANCE)
414 {
415  CGSolver cg;
416  cg.SetPrintLevel(print_iter);
417  cg.SetMaxIter(max_num_iter);
418  cg.SetRelTol(sqrt(RTOLERANCE));
419  cg.SetAbsTol(sqrt(ATOLERANCE));
420  cg.SetOperator(A);
421  cg.Mult(b, x);
422 }
423 
424 void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x,
425  int print_iter, int max_num_iter,
426  double RTOLERANCE, double ATOLERANCE)
427 {
428  CGSolver pcg;
429  pcg.SetPrintLevel(print_iter);
430  pcg.SetMaxIter(max_num_iter);
431  pcg.SetRelTol(sqrt(RTOLERANCE));
432  pcg.SetAbsTol(sqrt(ATOLERANCE));
433  pcg.SetOperator(A);
434  pcg.SetPreconditioner(B);
435  pcg.Mult(b, x);
436 }
437 
438 
439 inline void GeneratePlaneRotation(double &dx, double &dy,
440  double &cs, double &sn)
441 {
442  if (dy == 0.0)
443  {
444  cs = 1.0;
445  sn = 0.0;
446  }
447  else if (fabs(dy) > fabs(dx))
448  {
449  double temp = dx / dy;
450  sn = 1.0 / sqrt( 1.0 + temp*temp );
451  cs = temp * sn;
452  }
453  else
454  {
455  double temp = dy / dx;
456  cs = 1.0 / sqrt( 1.0 + temp*temp );
457  sn = temp * cs;
458  }
459 }
460 
461 inline void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
462 {
463  double temp = cs * dx + sn * dy;
464  dy = -sn * dx + cs * dy;
465  dx = temp;
466 }
467 
468 inline void Update(Vector &x, int k, DenseMatrix &h, Vector &s,
469  Array<Vector*> &v)
470 {
471  Vector y(s);
472 
473  // Backsolve:
474  for (int i = k; i >= 0; i--)
475  {
476  y(i) /= h(i,i);
477  for (int j = i - 1; j >= 0; j--)
478  y(j) -= h(j,i) * y(i);
479  }
480 
481  for (int j = 0; j <= k; j++)
482  x.Add(y(j), *v[j]);
483 }
484 
485 void GMRESSolver::Mult(const Vector &b, Vector &x) const
486 {
487  // Generalized Minimum Residual method following the algorithm
488  // on p. 20 of the SIAM Templates book.
489 
490  int n = width;
491 
492  DenseMatrix H(m+1, m);
493  Vector s(m+1), cs(m+1), sn(m+1);
494  Vector r(n), w(n);
495 
496  double resid;
497  int i, j, k;
498 
499  if (iterative_mode)
500  oper->Mult(x, r);
501  else
502  x = 0.0;
503 
504  if (prec)
505  {
506  if (iterative_mode)
507  {
508  subtract(b, r, w);
509  prec->Mult(w, r); // r = M (b - A x)
510  }
511  else
512  prec->Mult(b, r);
513  }
514  else
515  {
516  if (iterative_mode)
517  subtract(b, r, r);
518  else
519  r = b;
520  }
521  double beta = Norm(r); // beta = ||r||
522 
523  final_norm = fmax(rel_tol*beta, abs_tol);
524 
525  if (beta <= final_norm)
526  {
527  final_norm = beta;
528  final_iter = 0;
529  converged = 1;
530  return;
531  }
532 
533  if (print_level >= 0)
534  cout << " Pass : " << setw(2) << 1
535  << " Iteration : " << setw(3) << 0
536  << " ||B r|| = " << beta << '\n';
537 
538  Array<Vector *> v(m+1);
539  v = NULL;
540 
541  for (j = 1; j <= max_iter; )
542  {
543  if (v[0] == NULL) v[0] = new Vector(n);
544  v[0]->Set(1.0/beta, r);
545  s = 0.0; s(0) = beta;
546 
547  for (i = 0; i < m && j <= max_iter; i++, j++)
548  {
549  if (prec)
550  {
551  oper->Mult(*v[i], r);
552  prec->Mult(r, w); // w = M A v[i]
553  }
554  else
555  {
556  oper->Mult(*v[i], w);
557  }
558 
559  for (k = 0; k <= i; k++)
560  {
561  H(k,i) = Dot(w, *v[k]); // H(k,i) = w * v[k]
562  w.Add(-H(k,i), *v[k]); // w -= H(k,i) * v[k]
563  }
564 
565  H(i+1,i) = Norm(w); // H(i+1,i) = ||w||
566  if (v[i+1] == NULL) v[i+1] = new Vector(n);
567  v[i+1]->Set(1.0/H(i+1,i), w); // v[i+1] = w / H(i+1,i)
568 
569  for (k = 0; k < i; k++)
570  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
571 
572  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
573  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
574  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
575 
576  resid = fabs(s(i+1));
577  if (print_level >= 0)
578  cout << " Pass : " << setw(2) << (j-1)/m+1
579  << " Iteration : " << setw(3) << j
580  << " ||B r|| = " << resid << '\n';
581 
582  if (resid <= final_norm)
583  {
584  Update(x, i, H, s, v);
585  final_norm = resid;
586  final_iter = j;
587  for (i = 0; i <= m; i++)
588  delete v[i];
589  converged = 1;
590  return;
591  }
592  }
593 
594  if (print_level >= 0 && j <= max_iter)
595  cout << "Restarting..." << '\n';
596 
597  Update(x, i-1, H, s, v);
598 
599  oper->Mult(x, r);
600  if (prec)
601  {
602  subtract(b, r, w);
603  prec->Mult(w, r); // r = M (b - A x)
604  }
605  else
606  {
607  subtract(b, r, r);
608  }
609  beta = Norm(r); // beta = ||r||
610  if (beta <= final_norm)
611  {
612  final_norm = beta;
613  final_iter = j;
614  for (i = 0; i <= m; i++)
615  delete v[i];
616  converged = 1;
617  return;
618  }
619  }
620 
621  final_norm = beta;
623  for (i = 0; i <= m; i++)
624  delete v[i];
625  converged = 0;
626 }
627 
628 void FGMRESSolver::Mult(const Vector &b, Vector &x) const
629 {
630  DenseMatrix H(m+1,m);
631  Vector s(m+1), cs(m+1), sn(m+1);
632  Vector av(b.Size());
633  Vector r(b.Size());
634 
635  int i, j, k;
636 
637 
638  if (iterative_mode)
639  {
640  oper->Mult(x, r);
641  subtract(b,r,r);
642  }
643  else
644  {
645  x = 0.;
646  r = b;
647  }
648  double beta = Norm(r); // beta = ||r||
649 
650  final_norm = fmax(rel_tol*beta, abs_tol);
651 
652  if (beta <= final_norm)
653  {
654  final_norm = beta;
655  final_iter = 0;
656  converged = 1;
657  return;
658  }
659 
660  if (print_level>=0)
661  cout << " Pass : " << setw(2) << 1
662  << " Iteration : " << setw(3) << 0
663  << " || r || = " << beta << endl;
664 
665  Array<Vector*> v(m+1);
666  Array<Vector*> z(m+1);
667  for (i= 0; i<=m; i++)
668  {
669  v[i] = NULL;
670  z[i] = NULL;
671  }
672 
673  j = 1;
674  while (j <= max_iter)
675  {
676  if (v[0] == NULL) v[0] = new Vector(b.Size());
677  (*v[0]) = 0.0;
678  v[0] -> Add (1.0/beta, r); // v[0] = r / ||r||
679  s = 0.0; s(0) = beta;
680 
681  for (i = 0; i < m && j <= max_iter; i++)
682  {
683 
684  if (z[i] == NULL) z[i] = new Vector(b.Size());
685  (*z[i]) = 0.0;
686 
687  prec->Mult(*v[i], *z[i]);
688  oper->Mult(*z[i], r);
689 
690  for (k = 0; k <= i; k++)
691  {
692  H(k,i) = Dot( r, *v[k]); // H(k,i) = r * v[k]
693  r.Add(-H(k,i), (*v[k])); // r -= H(k,i) * v[k]
694  }
695 
696  H(i+1,i) = Norm(r); // H(i+1,i) = ||r||
697  if (v[i+1] == NULL) v[i+1] = new Vector(b.Size());
698  (*v[i+1]) = 0.0;
699  v[i+1] -> Add (1.0/H(i+1,i), r); // v[i+1] = r / H(i+1,i)
700 
701  for (k = 0; k < i; k++)
702  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
703 
704  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
705  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
706  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
707 
708  double resid = fabs(s(i+1));
709  if (print_level >= 0)
710  cout << " Pass : " << setw(2) << j
711  << " Iteration : " << setw(3) << i+1
712  << " || r || = " << resid << endl;
713 
714  if ( resid <= final_norm) {
715  Update(x, i, H, s, z);
716  final_norm = resid;
717  final_iter = (j-1)*m + i;
718  converged = 1;
719  for (i= 0; i<=m; i++){
720  if (v[i]) delete v[i];
721  if (z[i]) delete z[i];
722  }
723  return;
724  }
725  }
726 
727  if (print_level>=0)
728  cout << "Restarting..." << endl;
729 
730  Update(x, i-1, H, s, z);
731 
732  oper->Mult(x, r);
733  subtract(b,r,r);
734  beta = Norm(r);
735  if ( beta <= final_norm)
736  {
737  final_norm = beta;
738  final_iter = j*m;
739  converged = 1;
740  for (i= 0; i<=m; i++){
741  if (v[i]) delete v[i];
742  if (z[i]) delete z[i];
743  }
744  return;
745  }
746 
747  j++;
748  }
749 
750  for (i = 0; i <= m; i++)
751  {
752  if (v[i]) delete v[i];
753  if (z[i]) delete z[i];
754  }
755  converged = 0;
756  return;
757 
758 }
759 
760 
761 int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M,
762  int &max_iter, int m, double &tol, double atol, int printit)
763 {
764  GMRESSolver gmres;
765  gmres.SetPrintLevel(printit);
766  gmres.SetMaxIter(max_iter);
767  gmres.SetKDim(m);
768  gmres.SetRelTol(sqrt(tol));
769  gmres.SetAbsTol(sqrt(atol));
770  gmres.SetOperator(A);
771  gmres.SetPreconditioner(M);
772  gmres.Mult(b, x);
773  max_iter = gmres.GetNumIterations();
774  tol = gmres.GetFinalNorm()*gmres.GetFinalNorm();
775  return gmres.GetConverged();
776 }
777 
778 void GMRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
779  int print_iter, int max_num_iter, int m, double rtol, double atol)
780 {
781  GMRES(A, x, b, B, max_num_iter, m, rtol, atol, print_iter);
782 }
783 
784 
786 {
787  p.SetSize(width);
788  phat.SetSize(width);
789  s.SetSize(width);
790  shat.SetSize(width);
791  t.SetSize(width);
792  v.SetSize(width);
793  r.SetSize(width);
795 }
796 
797 void BiCGSTABSolver::Mult(const Vector &b, Vector &x) const
798 {
799  // BiConjugate Gradient Stabilized method following the algorithm
800  // on p. 27 of the SIAM Templates book.
801 
802  int i;
803  double resid, tol_goal;
804  double rho_1, rho_2=1.0, alpha=1.0, beta, omega=1.0;
805 
806  if (iterative_mode)
807  {
808  oper->Mult(x, r);
809  subtract(b, r, r); // r = b - A x
810  }
811  else
812  {
813  x = 0.0;
814  r = b;
815  }
816  rtilde = r;
817 
818  resid = Norm(r);
819  if (print_level >= 0)
820  cout << " Iteration : " << setw(3) << 0
821  << " ||r|| = " << resid << '\n';
822 
823  tol_goal = fmax(resid*rel_tol, abs_tol);
824 
825  if (resid <= tol_goal)
826  {
827  final_norm = resid;
828  final_iter = 0;
829  converged = 1;
830  return;
831  }
832 
833  for (i = 1; i <= max_iter; i++)
834  {
835  rho_1 = Dot(rtilde, r);
836  if (rho_1 == 0)
837  {
838  if (print_level >= 0)
839  cout << " Iteration : " << setw(3) << i
840  << " ||r|| = " << resid << '\n';
841  final_norm = resid;
842  final_iter = i;
843  converged = 0;
844  return;
845  }
846  if (i == 1)
847  p = r;
848  else
849  {
850  beta = (rho_1/rho_2) * (alpha/omega);
851  add(p, -omega, v, p); // p = p - omega * v
852  add(r, beta, p, p); // p = r + beta * p
853  }
854  if (prec)
855  prec->Mult(p, phat); // phat = M^{-1} * p
856  else
857  phat = p;
858  oper->Mult(phat, v); // v = A * phat
859  alpha = rho_1 / Dot(rtilde, v);
860  add(r, -alpha, v, s); // s = r - alpha * v
861  resid = Norm(s);
862  if (resid < tol_goal)
863  {
864  x.Add(alpha, phat); // x = x + alpha * phat
865  if (print_level >= 0)
866  cout << " Iteration : " << setw(3) << i
867  << " ||s|| = " << resid << '\n';
868  final_norm = resid;
869  final_iter = i;
870  converged = 1;
871  return;
872  }
873  if (print_level >= 0)
874  cout << " Iteration : " << setw(3) << i
875  << " ||s|| = " << resid;
876  if (prec)
877  prec->Mult(s, shat); // shat = M^{-1} * s
878  else
879  shat = s;
880  oper->Mult(shat, t); // t = A * shat
881  omega = Dot(t, s) / Dot(t, t);
882  x.Add(alpha, phat); // x += alpha * phat
883  x.Add(omega, shat); // x += omega * shat
884  add(s, -omega, t, r); // r = s - omega * t
885 
886  rho_2 = rho_1;
887  resid = Norm(r);
888  if (print_level >= 0)
889  cout << " ||r|| = " << resid << '\n';
890  if (resid < tol_goal)
891  {
892  final_norm = resid;
893  final_iter = i;
894  converged = 1;
895  return;
896  }
897  if (omega == 0)
898  {
899  final_norm = resid;
900  final_iter = i;
901  converged = 0;
902  return;
903  }
904  }
905 
906  final_norm = resid;
908  converged = 0;
909 }
910 
911 int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M,
912  int &max_iter, double &tol, double atol, int printit)
913 {
914  BiCGSTABSolver bicgstab;
915  bicgstab.SetPrintLevel(printit);
916  bicgstab.SetMaxIter(max_iter);
917  bicgstab.SetRelTol(sqrt(tol));
918  bicgstab.SetAbsTol(sqrt(atol));
919  bicgstab.SetOperator(A);
920  bicgstab.SetPreconditioner(M);
921  bicgstab.Mult(b, x);
922  max_iter = bicgstab.GetNumIterations();
923  tol = bicgstab.GetFinalNorm()*bicgstab.GetFinalNorm();
924  return bicgstab.GetConverged();
925 }
926 
927 void BiCGSTAB(const Operator &A, Solver &B, const Vector &b, Vector &x,
928  int print_iter, int max_num_iter, double rtol, double atol)
929 {
930  BiCGSTAB(A, x, b, B, max_num_iter, rtol, atol, print_iter);
931 }
932 
933 
935 {
937  v0.SetSize(width);
938  v1.SetSize(width);
939  w0.SetSize(width);
940  w1.SetSize(width);
941  q.SetSize(width);
942  if (prec)
943  u1.SetSize(width);
944 }
945 
946 void MINRESSolver::Mult(const Vector &b, Vector &x) const
947 {
948  // Based on the MINRES algorithm on p. 86, Fig. 6.9 in
949  // "Iterative Krylov Methods for Large Linear Systems",
950  // by Henk A. van der Vorst, 2003.
951  // Extended to support an SPD preconditioner.
952 
953  int it;
954  double beta, eta, gamma0, gamma1, sigma0, sigma1;
955  double alpha, delta, rho1, rho2, rho3, norm_goal;
956  Vector *z = (prec) ? &u1 : &v1;
957 
958  converged = 1;
959 
960  if (!iterative_mode)
961  {
962  v1 = b;
963  x = 0.;
964  }
965  else
966  {
967  oper->Mult(x, v1);
968  subtract(b, v1, v1);
969  }
970 
971  if (prec)
972  prec->Mult(v1, u1);
973  eta = beta = sqrt(Dot(*z, v1));
974  gamma0 = gamma1 = 1.;
975  sigma0 = sigma1 = 0.;
976 
977  norm_goal = fmax(rel_tol*eta, abs_tol);
978 
979  if (print_level == 1 || print_level == 3)
980  cout << "MINRES: iteration " << setw(3) << 0 << ": ||r||_B = "
981  << eta << '\n';
982 
983  if (eta <= norm_goal)
984  {
985  it = 0;
986  goto loop_end;
987  }
988 
989  for (it = 1; it <= max_iter; it++)
990  {
991  v1 /= beta;
992  if (prec)
993  u1 /= beta;
994  oper->Mult(*z, q);
995  alpha = Dot(*z, q);
996  if (it > 1) // (v0 == 0) for (it == 1)
997  q.Add(-beta, v0);
998  add(q, -alpha, v1, v0);
999 
1000  delta = gamma1*alpha - gamma0*sigma1*beta;
1001  rho3 = sigma0*beta;
1002  rho2 = sigma1*alpha + gamma0*gamma1*beta;
1003  if (!prec)
1004  {
1005  beta = Norm(v0);
1006  }
1007  else
1008  {
1009  prec->Mult(v0, q);
1010  beta = sqrt(Dot(v0, q));
1011  }
1012  rho1 = hypot(delta, beta);
1013 
1014  if (it == 1)
1015  w0.Set(1./rho1, *z); // (w0 == 0) and (w1 == 0)
1016  else if (it == 2)
1017  add(1./rho1, *z, -rho2/rho1, w1, w0); // (w0 == 0)
1018  else
1019  {
1020  add(-rho3/rho1, w0, -rho2/rho1, w1, w0);
1021  w0.Add(1./rho1, *z);
1022  }
1023 
1024  gamma0 = gamma1;
1025  gamma1 = delta/rho1;
1026 
1027  x.Add(gamma1*eta, w0);
1028 
1029  sigma0 = sigma1;
1030  sigma1 = beta/rho1;
1031 
1032  eta = -sigma1*eta;
1033 
1034  if (print_level == 1)
1035  cout << "MINRES: iteration " << setw(3) << it << ": ||r||_B = "
1036  << fabs(eta) << '\n';
1037 
1038  if (fabs(eta) <= norm_goal)
1039  goto loop_end;
1040 
1041  if (prec)
1042  swap(&u1, &q);
1043  swap(&v0, &v1);
1044  swap(&w0, &w1);
1045  }
1046  converged = 0;
1047  it--;
1048 
1049 loop_end:
1050  final_iter = it;
1051  final_norm = fabs(eta);
1052 
1053  if (print_level == 2)
1054  cout << "MINRES: number of iterations: " << it << '\n';
1055  else if (print_level == 3)
1056  cout << "MINRES: iteration " << setw(3) << it << ": ||r||_B = "
1057  << fabs(eta) << '\n';
1058 #if 0
1059  if (print_level >= 1)
1060  {
1061  oper->Mult(x, v1);
1062  subtract(b, v1, v1);
1063  if (prec)
1064  prec->Mult(v1, u1);
1065  eta = sqrt(Dot(*z, v1));
1066  cout << "MINRES: iteration " << setw(3) << it << ": ||r||_B = "
1067  << eta << " (re-computed)" << '\n';
1068  }
1069 #endif
1070  if (!converged && print_level >= 0)
1071  cerr << "MINRES: No convergence!" << '\n';
1072 }
1073 
1074 void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it,
1075  int max_it, double rtol, double atol)
1076 {
1077  MINRESSolver minres;
1078  minres.SetPrintLevel(print_it);
1079  minres.SetMaxIter(max_it);
1080  minres.SetRelTol(sqrt(rtol));
1081  minres.SetAbsTol(sqrt(atol));
1082  minres.SetOperator(A);
1083  minres.Mult(b, x);
1084 }
1085 
1086 void MINRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
1087  int print_it, int max_it, double rtol, double atol)
1088 {
1089  MINRESSolver minres;
1090  minres.SetPrintLevel(print_it);
1091  minres.SetMaxIter(max_it);
1092  minres.SetRelTol(sqrt(rtol));
1093  minres.SetAbsTol(sqrt(atol));
1094  minres.SetOperator(A);
1095  minres.SetPreconditioner(B);
1096  minres.Mult(b, x);
1097 }
1098 
1099 
1101 {
1102  oper = &op;
1103  height = op.Height();
1104  width = op.Width();
1105  MFEM_ASSERT(height == width, "square Operator is required.");
1106 
1107  r.SetSize(width);
1108  c.SetSize(width);
1109 }
1110 
1111 void NewtonSolver::Mult(const Vector &b, Vector &x) const
1112 {
1113  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
1114  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
1115 
1116  int it;
1117  double norm, norm_goal;
1118  bool have_b = (b.Size() == Height());
1119 
1120  if (!iterative_mode)
1121  x = 0.0;
1122 
1123  oper->Mult(x, r);
1124  if (have_b)
1125  r -= b;
1126 
1127  norm = Norm(r);
1128  norm_goal = fmax(rel_tol*norm, abs_tol);
1129 
1130  prec->iterative_mode = false;
1131 
1132  // x_{i+1} = x_i - [DF(x_i)]^{-1} [F(x_i)-b]
1133  for (it = 0; true; it++)
1134  {
1135  if (print_level >= 0)
1136  cout << "Newton iteration " << setw(2) << it
1137  << " : ||r|| = " << norm << '\n';
1138 
1139  if (norm <= norm_goal)
1140  {
1141  converged = 1;
1142  break;
1143  }
1144 
1145  if (it >= max_iter)
1146  {
1147  converged = 0;
1148  break;
1149  }
1150 
1152 
1153  prec->Mult(r, c); // c = [DF(x_i)]^{-1} [F(x_i)-b]
1154 
1155  x -= c;
1156 
1157  oper->Mult(x, r);
1158  if (have_b)
1159  r -= b;
1160  norm = Norm(r);
1161  }
1162 
1163  final_iter = it;
1164  final_norm = norm;
1165 }
1166 
1167 
1168 int aGMRES(const Operator &A, Vector &x, const Vector &b,
1169  const Operator &M, int &max_iter,
1170  int m_max, int m_min, int m_step, double cf,
1171  double &tol, double &atol, int printit)
1172 {
1173  int n = A.Width();
1174 
1175  int m = m_max;
1176 
1177  DenseMatrix H(m+1,m);
1178  Vector s(m+1), cs(m+1), sn(m+1);
1179  Vector w(n), av(n);
1180 
1181  double r1, resid;
1182  int i, j, k;
1183 
1184  M.Mult(b,w);
1185  double normb = w.Norml2(); // normb = ||M b||
1186  if (normb == 0.0)
1187  normb = 1;
1188 
1189  Vector r(n);
1190  A.Mult(x, r);
1191  subtract(b,r,w);
1192  M.Mult(w, r); // r = M (b - A x)
1193  double beta = r.Norml2(); // beta = ||r||
1194 
1195  resid = beta / normb;
1196 
1197  if (resid * resid <= tol) {
1198  tol = resid * resid;
1199  max_iter = 0;
1200  return 0;
1201  }
1202 
1203  if (printit)
1204  cout << " Pass : " << setw(2) << 1
1205  << " Iteration : " << setw(3) << 0
1206  << " (r, r) = " << beta*beta << '\n';
1207 
1208  tol *= (normb*normb);
1209  tol = (atol > tol) ? atol : tol;
1210 
1211  m = m_max;
1212  Array<Vector *> v(m+1);
1213  for (i= 0; i<=m; i++) {
1214  v[i] = new Vector(n);
1215  (*v[i]) = 0.0;
1216  }
1217 
1218  j = 1;
1219  while (j <= max_iter) {
1220  (*v[0]) = 0.0;
1221  v[0] -> Add (1.0/beta, r); // v[0] = r / ||r||
1222  s = 0.0; s(0) = beta;
1223 
1224  r1 = beta;
1225 
1226  for (i = 0; i < m && j <= max_iter; i++) {
1227  A.Mult((*v[i]),av);
1228  M.Mult(av,w); // w = M A v[i]
1229 
1230  for (k = 0; k <= i; k++) {
1231  H(k,i) = w * (*v[k]); // H(k,i) = w * v[k]
1232  w.Add(-H(k,i), (*v[k])); // w -= H(k,i) * v[k]
1233  }
1234 
1235  H(i+1,i) = w.Norml2(); // H(i+1,i) = ||w||
1236  (*v[i+1]) = 0.0;
1237  v[i+1] -> Add (1.0/H(i+1,i), w); // v[i+1] = w / H(i+1,i)
1238 
1239  for (k = 0; k < i; k++)
1240  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
1241 
1242  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1243  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1244  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
1245 
1246  resid = fabs(s(i+1));
1247  if (printit)
1248  cout << " Pass : " << setw(2) << j
1249  << " Iteration : " << setw(3) << i+1
1250  << " (r, r) = " << resid*resid << '\n';
1251 
1252  if ( resid*resid < tol) {
1253  Update(x, i, H, s, v);
1254  tol = resid * resid;
1255  max_iter = j;
1256  for (i= 0; i<=m; i++)
1257  delete v[i];
1258  return 0;
1259  }
1260  }
1261 
1262  if (printit)
1263  cout << "Restarting..." << '\n';
1264 
1265  Update(x, i-1, H, s, v);
1266 
1267  A.Mult(x, r);
1268  subtract(b,r,w);
1269  M.Mult(w, r); // r = M (b - A x)
1270  beta = r.Norml2(); // beta = ||r||
1271  if ( resid*resid < tol) {
1272  tol = resid * resid;
1273  max_iter = j;
1274  for (i= 0; i<=m; i++)
1275  delete v[i];
1276  return 0;
1277  }
1278 
1279  if (beta/r1 > cf)
1280  {
1281  if (m - m_step >= m_min)
1282  m -= m_step;
1283  else
1284  m = m_max;
1285  }
1286 
1287  j++;
1288  }
1289 
1290  tol = resid * resid;
1291  for (i= 0; i<=m; i++)
1292  delete v[i];
1293  return 1;
1294 }
1295 
1296 
1297 void SLBQPOptimizer::SetBounds(const Vector &_lo, const Vector &_hi)
1298 {
1299  lo.SetDataAndSize(_lo.GetData(), _lo.Size());
1300  hi.SetDataAndSize(_hi.GetData(), _hi.Size());
1301 }
1302 
1304 {
1305  w.SetDataAndSize(_w.GetData(), _w.Size());
1306  a = _a;
1307 }
1308 
1310 {
1311  mfem_error("SLBQPOptimizer::SetPreconditioner() : "
1312  "not meaningful for this solver");
1313 }
1314 
1316 {
1317  mfem_error("SLBQPOptimizer::SetOperator() : "
1318  "not meaningful for this solver");
1319 }
1320 
1321 inline void SLBQPOptimizer::print_iteration(int it, double r, double l) const
1322 {
1323  if (print_level > 1)
1324  cout << "SLBQP iteration " << it << ": residual = " << r
1325  << ", lambda = " << l << '\n';
1326 }
1327 
1328 void SLBQPOptimizer::Mult(const Vector& xt, Vector& x) const
1329 {
1330  // Based on code provided by Denis Ridzal, dridzal@sandia.gov.
1331  // Algorithm adapted from Dai and Fletcher, "New Algorithms for
1332  // Singly Linearly Constrained Quadratic Programs Subject to Lower
1333  // and Upper Bounds", Numerical Analysis Report NA/216, 2003.
1334 
1335  // Set some algorithm-specific constants and temporaries.
1336  int nclip = 0;
1337  double l = 0;
1338  double llow = 0;
1339  double lupp = 0;
1340  double lnew = 0;
1341  double dl = 2;
1342  double r = 0;
1343  double rlow = 0;
1344  double rupp = 0;
1345  double s = 0;
1346 
1347  const double smin = 0.1;
1348 
1349  const double tol = max(abs_tol, rel_tol*a);
1350 
1351  // *** Start bracketing phase of SLBQP ***
1352  if (print_level > 1)
1353  cout << "SLBQP bracketing phase" << '\n';
1354 
1355  // Solve QP with fixed Lagrange multiplier
1356  r = solve(l,xt,x,nclip);
1357  print_iteration(nclip, r, l);
1358 
1359 
1360  // If x=xt was already within bounds and satisfies the linear
1361  // constraint, then we already have the solution.
1362  if (fabs(r) <= tol)
1363  {
1364  converged = true;
1365  goto slbqp_done;
1366  }
1367 
1368  if (r < 0)
1369  {
1370  llow = l; rlow = r; l = l + dl;
1371 
1372  // Solve QP with fixed Lagrange multiplier
1373  r = solve(l,xt,x,nclip);
1374  print_iteration(nclip, r, l);
1375 
1376  while ((r < 0) && (nclip < max_iter))
1377  {
1378  llow = l;
1379  s = rlow/r - 1.0;
1380  if (s < smin) { s = smin; }
1381  dl = dl + dl/s;
1382  l = l + dl;
1383 
1384  // Solve QP with fixed Lagrange multiplier
1385  r = solve(l,xt,x,nclip);
1386  print_iteration(nclip, r, l);
1387  }
1388 
1389  lupp = l; rupp = r;
1390  }
1391  else
1392  {
1393  lupp = l; rupp = r; l = l - dl;
1394 
1395  // Solve QP with fixed Lagrange multiplier
1396  r = solve(l,xt,x,nclip);
1397  print_iteration(nclip, r, l);
1398 
1399  while ((r > 0) && (nclip < max_iter))
1400  {
1401  lupp = l;
1402  s = rupp/r - 1.0;
1403  if (s < smin) { s = smin; }
1404  dl = dl + dl/s;
1405  l = l - dl;
1406 
1407  // Solve QP with fixed Lagrange multiplier
1408  r = solve(l,xt,x,nclip);
1409  print_iteration(nclip, r, l);
1410  }
1411 
1412  llow = l; rlow = r;
1413  }
1414 
1415  // *** Stop bracketing phase of SLBQP ***
1416 
1417 
1418  // *** Start secant phase of SLBQP ***
1419  if (print_level > 1)
1420  cout << "SLBQP secant phase" << '\n';
1421 
1422  s = 1.0 - rlow/rupp; dl = dl/s; l = lupp - dl;
1423 
1424  // Solve QP with fixed Lagrange multiplier
1425  r = solve(l,xt,x,nclip);
1426  print_iteration(nclip, r, l);
1427 
1428  while ( (fabs(r) > tol) && (nclip < max_iter) )
1429  {
1430  if (r > 0)
1431  {
1432  if (s <= 2.0)
1433  {
1434  lupp = l; rupp = r; s = 1.0 - rlow/rupp;
1435  dl = (lupp - llow)/s; l = lupp - dl;
1436  }
1437  else
1438  {
1439  s = rupp/r - 1.0;
1440  if (s < smin) { s = smin; }
1441  dl = (lupp - l)/s;
1442  lnew = 0.75*llow + 0.25*l;
1443  if (lnew < l-dl) { lnew = l-dl; }
1444  lupp = l; rupp = r; l = lnew;
1445  s = (lupp - llow)/(lupp - l);
1446  }
1447 
1448  }
1449  else
1450  {
1451  if (s >= 2.0)
1452  {
1453  llow = l; rlow = r; s = 1.0 - rlow/rupp;
1454  dl = (lupp - llow)/s; l = lupp - dl;
1455  }
1456  else
1457  {
1458  s = rlow/r - 1.0;
1459  if (s < smin) { s = smin; }
1460  dl = (l - llow)/s;
1461  lnew = 0.75*lupp + 0.25*l;
1462  if (lnew < l+dl) { lnew = l+dl; }
1463  llow = l; rlow = r; l = lnew;
1464  s = (lupp - llow)/(lupp - l);
1465  }
1466  }
1467 
1468  // Solve QP with fixed Lagrange multiplier
1469  r = solve(l,xt,x,nclip);
1470  print_iteration(nclip, r, l);
1471  }
1472 
1473  // *** Stop secant phase of SLBQP ***
1474 
1475  converged = (fabs(r) <= tol);
1476  if (!converged && print_level >= 0)
1477  cerr << "SLBQP not converged!" << '\n';
1478 
1479 slbqp_done:
1480 
1481  final_iter = nclip;
1482  final_norm = r;
1483 
1484  if (print_level == 1 || (!converged && print_level >= 0))
1485  {
1486  cout << "SLBQP iterations = " << nclip << '\n';
1487  cout << "SLBQP lamba = " << l << '\n';
1488  cout << "SLBQP residual = " << r << '\n';
1489  }
1490 }
1491 
1492 #ifdef MFEM_USE_SUITESPARSE
1493 
1495 {
1496  mat = NULL;
1497  Numeric = NULL;
1498  AI = AJ = NULL;
1499  if (!use_long_ints)
1500  umfpack_di_defaults(Control);
1501  else
1502  umfpack_dl_defaults(Control);
1503 }
1504 
1506 {
1507  int *Ap, *Ai;
1508  void *Symbolic;
1509  double *Ax;
1510 
1511  if (Numeric)
1512  {
1513  if (!use_long_ints)
1514  {
1515  umfpack_di_free_numeric(&Numeric);
1516  }
1517  else
1518  {
1519  umfpack_dl_free_numeric(&Numeric);
1520  }
1521  }
1522 
1523  mat = const_cast<SparseMatrix *>(dynamic_cast<const SparseMatrix *>(&op));
1524  if (mat == NULL)
1525  mfem_error("UMFPackSolver::SetOperator : not a SparseMatrix!");
1526 
1527  // UMFPack requires that the column-indices in mat corresponding to each
1528  // row be sorted.
1529  // Generally, this will modify the ordering of the entries of mat.
1531 
1532  height = mat->Height();
1533  width = mat->Width();
1534  Ap = mat->GetI();
1535  Ai = mat->GetJ();
1536  Ax = mat->GetData();
1537 
1538  if (!use_long_ints)
1539  {
1540  int status = umfpack_di_symbolic(width, width, Ap, Ai, Ax, &Symbolic,
1541  Control, Info);
1542  if (status < 0)
1543  {
1544  umfpack_di_report_info(Control, Info);
1545  umfpack_di_report_status(Control, status);
1546  mfem_error("UMFPackSolver::SetOperator :"
1547  " umfpack_di_symbolic() failed!");
1548  }
1549 
1550  status = umfpack_di_numeric(Ap, Ai, Ax, Symbolic, &Numeric,
1551  Control, Info);
1552  if (status < 0)
1553  {
1554  umfpack_di_report_info(Control, Info);
1555  umfpack_di_report_status(Control, status);
1556  mfem_error("UMFPackSolver::SetOperator :"
1557  " umfpack_di_numeric() failed!");
1558  }
1559  umfpack_di_free_symbolic(&Symbolic);
1560  }
1561  else
1562  {
1563  SuiteSparse_long status;
1564 
1565  delete [] AJ;
1566  delete [] AI;
1567  AI = new SuiteSparse_long[width + 1];
1568  AJ = new SuiteSparse_long[Ap[width]];
1569  for (int i = 0; i <= width; i++)
1570  AI[i] = (SuiteSparse_long)(Ap[i]);
1571  for (int i = 0; i <= Ap[width]; i++)
1572  AJ[i] = (SuiteSparse_long)(Ai[i]);
1573 
1574  status = umfpack_dl_symbolic(width, width, AI, AJ, Ax, &Symbolic,
1575  Control, Info);
1576  if (status < 0)
1577  {
1578  umfpack_dl_report_info(Control, Info);
1579  umfpack_dl_report_status(Control, status);
1580  mfem_error("UMFPackSolver::SetOperator :"
1581  " umfpack_dl_symbolic() failed!");
1582  }
1583 
1584  status = umfpack_dl_numeric(AI, AJ, Ax, Symbolic, &Numeric,
1585  Control, Info);
1586  if (status < 0)
1587  {
1588  umfpack_dl_report_info(Control, Info);
1589  umfpack_dl_report_status(Control, status);
1590  mfem_error("UMFPackSolver::SetOperator :"
1591  " umfpack_dl_numeric() failed!");
1592  }
1593  umfpack_dl_free_symbolic(&Symbolic);
1594  }
1595 }
1596 
1597 void UMFPackSolver::Mult(const Vector &b, Vector &x) const
1598 {
1599  if (mat == NULL)
1600  mfem_error("UMFPackSolver::Mult : matrix is not set!"
1601  " Call SetOperator first!");
1602 
1603  if (!use_long_ints)
1604  {
1605  int status =
1606  umfpack_di_solve(UMFPACK_At, mat->GetI(), mat->GetJ(),
1607  mat->GetData(), x, b, Numeric, Control, Info);
1608  umfpack_di_report_info(Control, Info);
1609  if (status < 0)
1610  {
1611  umfpack_di_report_status(Control, status);
1612  mfem_error("UMFPackSolver::Mult : umfpack_di_solve() failed!");
1613  }
1614  }
1615  else
1616  {
1617  SuiteSparse_long status =
1618  umfpack_dl_solve(UMFPACK_At, AI, AJ, mat->GetData(), x, b,
1619  Numeric, Control, Info);
1620  umfpack_dl_report_info(Control, Info);
1621  if (status < 0)
1622  {
1623  umfpack_dl_report_status(Control, status);
1624  mfem_error("UMFPackSolver::Mult : umfpack_dl_solve() failed!");
1625  }
1626  }
1627 }
1628 
1630 {
1631  if (mat == NULL)
1632  mfem_error("UMFPackSolver::MultTranspose : matrix is not set!"
1633  " Call SetOperator first!");
1634 
1635  if (!use_long_ints)
1636  {
1637  int status =
1638  umfpack_di_solve(UMFPACK_A, mat->GetI(), mat->GetJ(),
1639  mat->GetData(), x, b, Numeric, Control, Info);
1640  umfpack_di_report_info(Control, Info);
1641  if (status < 0)
1642  {
1643  umfpack_di_report_status(Control, status);
1644  mfem_error("UMFPackSolver::MultTranspose :"
1645  " umfpack_di_solve() failed!");
1646  }
1647  }
1648  else
1649  {
1650  SuiteSparse_long status =
1651  umfpack_dl_solve(UMFPACK_A, AI, AJ, mat->GetData(), x, b,
1652  Numeric, Control, Info);
1653  umfpack_dl_report_info(Control, Info);
1654  if (status < 0)
1655  {
1656  umfpack_dl_report_status(Control, status);
1657  mfem_error("UMFPackSolver::MultTranspose :"
1658  " umfpack_dl_solve() failed!");
1659  }
1660  }
1661 }
1662 
1664 {
1665  delete [] AJ;
1666  delete [] AI;
1667  if (Numeric)
1668  {
1669  if (!use_long_ints)
1670  {
1671  umfpack_di_free_numeric(&Numeric);
1672  }
1673  else
1674  {
1675  umfpack_dl_free_numeric(&Numeric);
1676  }
1677  }
1678 }
1679 
1680 #endif // MFEM_USE_SUITESPARSE
1681 
1682 }
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: solvers.cpp:461
Conjugate gradient method.
Definition: solvers.hpp:109
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator.
Definition: solvers.cpp:1629
double Info[UMFPACK_INFO]
Definition: solvers.hpp:341
SparseMatrix * mat
Definition: solvers.hpp:333
const Operator * oper
Definition: solvers.hpp:39
SuiteSparse_long * AI
Definition: solvers.hpp:335
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:248
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x.
Definition: operator.hpp:54
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:532
Data type dense matrix.
Definition: densemat.hpp:22
int Size() const
Returns the size of the vector.
Definition: vector.hpp:76
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application.
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
Definition: operator.hpp:106
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:285
MINRES method.
Definition: solvers.hpp:218
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:112
double * GetData() const
Definition: vector.hpp:80
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:198
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:227
double Norm(const Vector &x) const
Definition: solvers.hpp:50
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:797
int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, double &tol, double atol, int printit)
BiCGSTAB method. (tolerances are squared)
Definition: solvers.cpp:911
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:485
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2430
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:946
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:628
void SetLinearConstraint(const Vector &_w, double _a)
Definition: solvers.cpp:1303
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:231
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:96
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Definition: solvers.hpp:296
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it, int max_it, double rtol, double atol)
MINRES method without preconditioner. (tolerances are squared)
Definition: solvers.cpp:1074
virtual ~UMFPackSolver()
Definition: solvers.cpp:1663
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:375
void SetMaxIter(int max_it)
Definition: solvers.hpp:61
SuiteSparse_long * AJ
Definition: solvers.hpp:335
void SetKDim(int dim)
Definition: solvers.hpp:153
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:411
void print_iteration(int it, double r, double l) const
Definition: solvers.cpp:1321
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:424
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:934
double Dot(const Vector &x, const Vector &y) const
Definition: solvers.cpp:50
double * GetData() const
Return element data.
Definition: sparsemat.hpp:101
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1111
int * GetI() const
Return the array I.
Definition: sparsemat.hpp:97
Stationary linear iteration: x &lt;- x + B (b - A x)
Definition: solvers.hpp:77
void SetAbsTol(double atol)
Definition: solvers.hpp:60
void SetRelTol(double rtol)
Definition: solvers.hpp:59
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:340
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:341
void swap(Vector *v1, Vector *v2)
Definition: vector.cpp:213
void mfem_error(const char *msg)
Definition: error.cpp:23
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1315
GMRES method.
Definition: solvers.hpp:141
void Update(Vector &x, int k, DenseMatrix &h, Vector &s, Array< Vector * > &v)
Definition: solvers.cpp:468
void SetDataAndSize(double *d, int s)
Definition: vector.hpp:64
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:182
void SetBounds(const Vector &_lo, const Vector &_hi)
Definition: solvers.cpp:1297
void UpdateVectors()
Definition: solvers.cpp:106
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:168
void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: solvers.cpp:439
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:201
virtual void Mult(const Vector &xt, Vector &x) const
Operator application.
Definition: solvers.cpp:1328
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:123
Vector data type.
Definition: vector.hpp:29
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:90
int aGMRES(const Operator &A, Vector &x, const Vector &b, const Operator &M, int &max_iter, int m_max, int m_min, int m_step, double cf, double &tol, double &atol, int printit)
Definition: solvers.cpp:1168
void UpdateVectors()
Definition: solvers.cpp:278
BiCGSTAB method.
Definition: solvers.hpp:187
Base class for solvers.
Definition: operator.hpp:102
double GetFinalNorm()
Definition: solvers.hpp:66
virtual void SetPreconditioner(Solver &pr)
These are not currently meaningful for this solver and will error out.
Definition: solvers.cpp:1309
Abstract operator.
Definition: operator.hpp:21
int * GetJ() const
Return the array J.
Definition: sparsemat.hpp:99
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1597
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:91
void SLI(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Stationary linear iteration. (tolerances are squared)
Definition: solvers.cpp:250
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
Definition: solvers.cpp:761
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1100
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1505