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