MFEM  v4.5.2
Finite element discretization library
solvers.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "linalg.hpp"
13 #include "../general/annotation.hpp"
14 #include "../general/forall.hpp"
15 #include "../general/globals.hpp"
16 #include "../fem/bilinearform.hpp"
17 #include <iostream>
18 #include <iomanip>
19 #include <algorithm>
20 #include <cmath>
21 #include <set>
22 
23 namespace mfem
24 {
25 
26 using namespace std;
27 
29  : Solver(0, true)
30 {
31  oper = NULL;
32  prec = NULL;
33  max_iter = 10;
34  rel_tol = abs_tol = 0.0;
35 #ifdef MFEM_USE_MPI
36  dot_prod_type = 0;
37 #endif
38 }
39 
40 #ifdef MFEM_USE_MPI
41 
43  : Solver(0, true)
44 {
45  oper = NULL;
46  prec = NULL;
47  max_iter = 10;
48  rel_tol = abs_tol = 0.0;
49  dot_prod_type = 1;
50  comm = comm_;
51 }
52 
53 #endif // MFEM_USE_MPI
54 
55 double IterativeSolver::Dot(const Vector &x, const Vector &y) const
56 {
57 #ifndef MFEM_USE_MPI
58  return (x * y);
59 #else
60  if (dot_prod_type == 0)
61  {
62  return (x * y);
63  }
64  else
65  {
66  return InnerProduct(comm, x, y);
67  }
68 #endif
69 }
70 
71 void IterativeSolver::SetPrintLevel(int print_lvl)
72 {
74  int print_level_ = print_lvl;
75 
76 #ifdef MFEM_USE_MPI
77  if (dot_prod_type != 0)
78  {
79  int rank;
80  MPI_Comm_rank(comm, &rank);
81  if (rank != 0) // Suppress output.
82  {
83  print_level_ = -1;
85  }
86  }
87 #endif
88 
89  print_level = print_level_;
90 }
91 
93 {
94  print_options = options;
95 
96  int derived_print_level = GuessLegacyPrintLevel(options);
97 
98 #ifdef MFEM_USE_MPI
99  if (dot_prod_type != 0)
100  {
101  int rank;
102  MPI_Comm_rank(comm, &rank);
103  if (rank != 0)
104  {
105  derived_print_level = -1;
107  }
108  }
109 #endif
110 
111  print_level = derived_print_level;
112 }
113 
115  int print_level_)
116 {
117 #ifdef MFEM_USE_MPI
118  int rank = 0;
119  if (comm != MPI_COMM_NULL)
120  {
121  MPI_Comm_rank(comm, &rank);
122  }
123 #endif
124 
125  switch (print_level_)
126  {
127  case -1:
128  return PrintLevel();
129  case 0:
130  return PrintLevel().Errors().Warnings();
131  case 1:
132  return PrintLevel().Errors().Warnings().Iterations();
133  case 2:
134  return PrintLevel().Errors().Warnings().Summary();
135  case 3:
136  return PrintLevel().Errors().Warnings().FirstAndLast();
137  default:
138 #ifdef MFEM_USE_MPI
139  if (rank == 0)
140 #endif
141  {
142  MFEM_WARNING("Unknown print level " << print_level_ <<
143  ". Defaulting to level 0.");
144  }
145  return PrintLevel().Errors().Warnings();
146  }
147 }
148 
150 {
151  if (print_options_.iterations)
152  {
153  return 1;
154  }
155  else if (print_options_.first_and_last)
156  {
157  return 3;
158  }
159  else if (print_options_.summary)
160  {
161  return 2;
162  }
163  else if (print_options_.errors && print_options_.warnings)
164  {
165  return 0;
166  }
167  else
168  {
169  return -1;
170  }
171 }
172 
174 {
175  prec = &pr;
176  prec->iterative_mode = false;
177 }
178 
180 {
181  oper = &op;
182  height = op.Height();
183  width = op.Width();
184  if (prec)
185  {
186  prec->SetOperator(*oper);
187  }
188 }
189 
190 void IterativeSolver::Monitor(int it, double norm, const Vector& r,
191  const Vector& x, bool final) const
192 {
193  if (monitor != nullptr)
194  {
195  monitor->MonitorResidual(it, norm, r, final);
196  monitor->MonitorSolution(it, norm, x, final);
197  }
198 }
199 
201  : damping(dmpng),
202  ess_tdof_list(nullptr),
203  oper(nullptr),
204  allow_updates(true)
205 { }
206 
208  const Array<int> &ess_tdofs,
209  const double dmpng)
210  :
211  Solver(a.FESpace()->GetTrueVSize()),
212  dinv(height),
213  damping(dmpng),
214  ess_tdof_list(&ess_tdofs),
215  residual(height),
216  allow_updates(false)
217 {
218  Vector &diag(residual);
219  a.AssembleDiagonal(diag);
220  // 'a' cannot be used for iterative_mode == true because its size may be
221  // different.
222  oper = nullptr;
223  Setup(diag);
224 }
225 
227  const Array<int> &ess_tdofs,
228  const double dmpng)
229  :
230  Solver(d.Size()),
231  dinv(height),
232  damping(dmpng),
233  ess_tdof_list(&ess_tdofs),
234  residual(height),
235  oper(NULL),
236  allow_updates(false)
237 {
238  Setup(d);
239 }
240 
242 {
243  if (!allow_updates)
244  {
245  // original behavior of this method
246  oper = &op; return;
247  }
248 
249  // Treat (Par)BilinearForm objects as a special case since their
250  // AssembleDiagonal method returns the true-dof diagonal whereas the form
251  // itself may act as an ldof operator. This is for compatibility with the
252  // constructor that takes a BilinearForm parameter.
253  const BilinearForm *blf = dynamic_cast<const BilinearForm *>(&op);
254  if (blf)
255  {
256  // 'a' cannot be used for iterative_mode == true because its size may be
257  // different.
258  oper = nullptr;
259  height = width = blf->FESpace()->GetTrueVSize();
260  }
261  else
262  {
263  oper = &op;
264  height = op.Height();
265  width = op.Width();
266  MFEM_ASSERT(height == width, "not a square matrix!");
267  // ess_tdof_list is only used with BilinearForm
268  ess_tdof_list = nullptr;
269  }
270  dinv.SetSize(height);
271  residual.SetSize(height);
272  Vector &diag(residual);
273  op.AssembleDiagonal(diag);
274  Setup(diag);
275 }
276 
278 {
279  residual.UseDevice(true);
280  const double delta = damping;
281  auto D = diag.Read();
282  auto DI = dinv.Write();
283  const bool use_abs_diag_ = use_abs_diag;
284  MFEM_FORALL(i, height,
285  {
286  if (D[i] == 0.0)
287  {
288  MFEM_ABORT_KERNEL("Zero diagonal entry in OperatorJacobiSmoother");
289  }
290  if (!use_abs_diag_) { DI[i] = delta / D[i]; }
291  else { DI[i] = delta / std::abs(D[i]); }
292  });
293  if (ess_tdof_list && ess_tdof_list->Size() > 0)
294  {
295  auto I = ess_tdof_list->Read();
296  MFEM_FORALL(i, ess_tdof_list->Size(), DI[I[i]] = delta; );
297  }
298 }
299 
301 {
302  // For empty MPI ranks, height may be 0:
303  // MFEM_VERIFY(Height() > 0, "The diagonal hasn't been computed.");
304  MFEM_ASSERT(x.Size() == Width(), "invalid input vector");
305  MFEM_ASSERT(y.Size() == Height(), "invalid output vector");
306 
307  if (iterative_mode)
308  {
309  MFEM_VERIFY(oper, "iterative_mode == true requires the forward operator");
310  oper->Mult(y, residual); // r = A y
311  subtract(x, residual, residual); // r = x - A y
312  }
313  else
314  {
315  residual = x;
316  y.UseDevice(true);
317  y = 0.0;
318  }
319  auto DI = dinv.Read();
320  auto R = residual.Read();
321  auto Y = y.ReadWrite();
322  MFEM_FORALL(i, height,
323  {
324  Y[i] += DI[i] * R[i];
325  });
326 }
327 
329  const Vector &d,
330  const Array<int>& ess_tdofs,
331  int order_, double max_eig_estimate_)
332  :
333  Solver(d.Size()),
334  order(order_),
335  max_eig_estimate(max_eig_estimate_),
336  N(d.Size()),
337  dinv(N),
338  diag(d),
339  coeffs(order),
340  ess_tdof_list(ess_tdofs),
341  residual(N),
342  oper(&oper_) { Setup(); }
343 
344 #ifdef MFEM_USE_MPI
346  const Vector &d,
347  const Array<int>& ess_tdofs,
348  int order_, MPI_Comm comm, int power_iterations, double power_tolerance)
349 #else
351  const Vector &d,
352  const Array<int>& ess_tdofs,
353  int order_, int power_iterations, double power_tolerance)
354 #endif
355  : Solver(d.Size()),
356  order(order_),
357  N(d.Size()),
358  dinv(N),
359  diag(d),
360  coeffs(order),
361  ess_tdof_list(ess_tdofs),
362  residual(N),
363  oper(&oper_)
364 {
365  OperatorJacobiSmoother invDiagOperator(diag, ess_tdofs, 1.0);
366  ProductOperator diagPrecond(&invDiagOperator, oper, false, false);
367 
368 #ifdef MFEM_USE_MPI
369  PowerMethod powerMethod(comm);
370 #else
371  PowerMethod powerMethod;
372 #endif
373  Vector ev(oper->Width());
374  max_eig_estimate = powerMethod.EstimateLargestEigenvalue(diagPrecond, ev,
375  power_iterations, power_tolerance);
376 
377  Setup();
378 }
379 
381  const Vector &d,
382  const Array<int>& ess_tdofs,
383  int order_, double max_eig_estimate_)
384  : OperatorChebyshevSmoother(*oper_, d, ess_tdofs, order_, max_eig_estimate_) { }
385 
386 #ifdef MFEM_USE_MPI
388  const Vector &d,
389  const Array<int>& ess_tdofs,
390  int order_, MPI_Comm comm, int power_iterations, double power_tolerance)
391  : OperatorChebyshevSmoother(*oper_, d, ess_tdofs, order_, comm,
392  power_iterations, power_tolerance) { }
393 #else
395  const Vector &d,
396  const Array<int>& ess_tdofs,
397  int order_, int power_iterations, double power_tolerance)
398  : OperatorChebyshevSmoother(*oper_, d, ess_tdofs, order_, power_iterations,
399  power_tolerance) { }
400 #endif
401 
403 {
404  // Invert diagonal
405  residual.UseDevice(true);
406  auto D = diag.Read();
407  auto X = dinv.Write();
408  MFEM_FORALL(i, N, X[i] = 1.0 / D[i]; );
409  auto I = ess_tdof_list.Read();
410  MFEM_FORALL(i, ess_tdof_list.Size(), X[I[i]] = 1.0; );
411 
412  // Set up Chebyshev coefficients
413  // For reference, see e.g., Parallel multigrid smoothing: polynomial versus
414  // Gauss-Seidel by Adams et al.
415  double upper_bound = 1.2 * max_eig_estimate;
416  double lower_bound = 0.3 * max_eig_estimate;
417  double theta = 0.5 * (upper_bound + lower_bound);
418  double delta = 0.5 * (upper_bound - lower_bound);
419 
420  switch (order-1)
421  {
422  case 0:
423  {
424  coeffs[0] = 1.0 / theta;
425  break;
426  }
427  case 1:
428  {
429  double tmp_0 = 1.0/(pow(delta, 2) - 2*pow(theta, 2));
430  coeffs[0] = -4*theta*tmp_0;
431  coeffs[1] = 2*tmp_0;
432  break;
433  }
434  case 2:
435  {
436  double tmp_0 = 3*pow(delta, 2);
437  double tmp_1 = pow(theta, 2);
438  double tmp_2 = 1.0/(-4*pow(theta, 3) + theta*tmp_0);
439  coeffs[0] = tmp_2*(tmp_0 - 12*tmp_1);
440  coeffs[1] = 12/(tmp_0 - 4*tmp_1);
441  coeffs[2] = -4*tmp_2;
442  break;
443  }
444  case 3:
445  {
446  double tmp_0 = pow(delta, 2);
447  double tmp_1 = pow(theta, 2);
448  double tmp_2 = 8*tmp_0;
449  double tmp_3 = 1.0/(pow(delta, 4) + 8*pow(theta, 4) - tmp_1*tmp_2);
450  coeffs[0] = tmp_3*(32*pow(theta, 3) - 16*theta*tmp_0);
451  coeffs[1] = tmp_3*(-48*tmp_1 + tmp_2);
452  coeffs[2] = 32*theta*tmp_3;
453  coeffs[3] = -8*tmp_3;
454  break;
455  }
456  case 4:
457  {
458  double tmp_0 = 5*pow(delta, 4);
459  double tmp_1 = pow(theta, 4);
460  double tmp_2 = pow(theta, 2);
461  double tmp_3 = pow(delta, 2);
462  double tmp_4 = 60*tmp_3;
463  double tmp_5 = 20*tmp_3;
464  double tmp_6 = 1.0/(16*pow(theta, 5) - pow(theta, 3)*tmp_5 + theta*tmp_0);
465  double tmp_7 = 160*tmp_2;
466  double tmp_8 = 1.0/(tmp_0 + 16*tmp_1 - tmp_2*tmp_5);
467  coeffs[0] = tmp_6*(tmp_0 + 80*tmp_1 - tmp_2*tmp_4);
468  coeffs[1] = tmp_8*(tmp_4 - tmp_7);
469  coeffs[2] = tmp_6*(-tmp_5 + tmp_7);
470  coeffs[3] = -80*tmp_8;
471  coeffs[4] = 16*tmp_6;
472  break;
473  }
474  default:
475  MFEM_ABORT("Chebyshev smoother not implemented for order = " << order);
476  }
477 }
478 
480 {
481  if (iterative_mode)
482  {
483  MFEM_ABORT("Chebyshev smoother not implemented for iterative mode");
484  }
485 
486  if (!oper)
487  {
488  MFEM_ABORT("Chebyshev smoother requires operator");
489  }
490 
491  residual = x;
492  helperVector.SetSize(x.Size());
493 
494  y.UseDevice(true);
495  y = 0.0;
496 
497  for (int k = 0; k < order; ++k)
498  {
499  // Apply
500  if (k > 0)
501  {
502  oper->Mult(residual, helperVector);
503  residual = helperVector;
504  }
505 
506  // Scale residual by inverse diagonal
507  const int n = N;
508  auto Dinv = dinv.Read();
509  auto R = residual.ReadWrite();
510  MFEM_FORALL(i, n, R[i] *= Dinv[i]; );
511 
512  // Add weighted contribution to y
513  auto Y = y.ReadWrite();
514  auto C = coeffs.Read();
515  MFEM_FORALL(i, n, Y[i] += C[k] * R[i]; );
516  }
517 }
518 
520 {
521  r.SetSize(width);
522  z.SetSize(width);
523 }
524 
525 void SLISolver::Mult(const Vector &b, Vector &x) const
526 {
527  int i;
528 
529  // Optimized preconditioned SLI with fixed number of iterations and given
530  // initial guess
531  if (rel_tol == 0.0 && iterative_mode && prec)
532  {
533  for (i = 0; i < max_iter; i++)
534  {
535  oper->Mult(x, r); // r = A x
536  subtract(b, r, r); // r = b - A x
537  prec->Mult(r, z); // z = B r
538  add(x, 1.0, z, x); // x = x + B (b - A x)
539  }
540  converged = true;
541  final_iter = i;
542  return;
543  }
544 
545  // Optimized preconditioned SLI with fixed number of iterations and zero
546  // initial guess
547  if (rel_tol == 0.0 && !iterative_mode && prec)
548  {
549  prec->Mult(b, x); // x = B b (initial guess 0)
550  for (i = 1; i < max_iter; i++)
551  {
552  oper->Mult(x, r); // r = A x
553  subtract(b, r, r); // r = b - A x
554  prec->Mult(r, z); // z = B r
555  add(x, 1.0, z, x); // x = x + B (b - A x)
556  }
557  converged = true;
558  final_iter = i;
559  return;
560  }
561 
562  // General version of SLI with a relative tolerance, optional preconditioner
563  // and optional initial guess
564  double r0, nom, nom0, nomold = 1, cf;
565 
566  if (iterative_mode)
567  {
568  oper->Mult(x, r);
569  subtract(b, r, r); // r = b - A x
570  }
571  else
572  {
573  r = b;
574  x = 0.0;
575  }
576 
577  if (prec)
578  {
579  prec->Mult(r, z); // z = B r
580  nom0 = nom = sqrt(Dot(z, z));
581  }
582  else
583  {
584  nom0 = nom = sqrt(Dot(r, r));
585  }
586  initial_norm = nom0;
587 
589  {
590  mfem::out << " Iteration : " << setw(3) << right << 0 << " ||Br|| = "
591  << nom << (print_options.first_and_last ? " ..." : "") << '\n';
592  }
593 
594  r0 = std::max(nom*rel_tol, abs_tol);
595  if (nom <= r0)
596  {
597  converged = true;
598  final_iter = 0;
599  final_norm = nom;
600  return;
601  }
602 
603  // start iteration
604  converged = false;
606  for (i = 1; true; )
607  {
608  if (prec) // x = x + B (b - A x)
609  {
610  add(x, 1.0, z, x);
611  }
612  else
613  {
614  add(x, 1.0, r, x);
615  }
616 
617  oper->Mult(x, r);
618  subtract(b, r, r); // r = b - A x
619 
620  if (prec)
621  {
622  prec->Mult(r, z); // z = B r
623  nom = sqrt(Dot(z, z));
624  }
625  else
626  {
627  nom = sqrt(Dot(r, r));
628  }
629 
630  cf = nom/nomold;
631  nomold = nom;
632 
633  bool done = false;
634  if (nom < r0)
635  {
636  converged = true;
637  final_iter = i;
638  done = true;
639  }
640 
641  if (++i > max_iter)
642  {
643  done = true;
644  }
645 
647  {
648  mfem::out << " Iteration : " << setw(3) << right << (i-1)
649  << " ||Br|| = " << setw(11) << left << nom
650  << "\tConv. rate: " << cf << '\n';
651  }
652 
653  if (done) { break; }
654  }
655 
657  {
658  const auto rf = pow (nom/nom0, 1.0/final_iter);
659  mfem::out << "SLI: Number of iterations: " << final_iter << '\n'
660  << "Conv. rate: " << cf << '\n'
661  << "Average reduction factor: "<< rf << '\n';
662  }
664  {
665  mfem::out << "SLI: No convergence!" << '\n';
666  }
667 
668  final_norm = nom;
669 }
670 
671 void SLI(const Operator &A, const Vector &b, Vector &x,
672  int print_iter, int max_num_iter,
673  double RTOLERANCE, double ATOLERANCE)
674 {
675  MFEM_PERF_FUNCTION;
676 
677  SLISolver sli;
678  sli.SetPrintLevel(print_iter);
679  sli.SetMaxIter(max_num_iter);
680  sli.SetRelTol(sqrt(RTOLERANCE));
681  sli.SetAbsTol(sqrt(ATOLERANCE));
682  sli.SetOperator(A);
683  sli.Mult(b, x);
684 }
685 
686 void SLI(const Operator &A, Solver &B, const Vector &b, Vector &x,
687  int print_iter, int max_num_iter,
688  double RTOLERANCE, double ATOLERANCE)
689 {
690  MFEM_PERF_FUNCTION;
691 
692  SLISolver sli;
693  sli.SetPrintLevel(print_iter);
694  sli.SetMaxIter(max_num_iter);
695  sli.SetRelTol(sqrt(RTOLERANCE));
696  sli.SetAbsTol(sqrt(ATOLERANCE));
697  sli.SetOperator(A);
698  sli.SetPreconditioner(B);
699  sli.Mult(b, x);
700 }
701 
702 
704 {
706 
707  r.SetSize(width, mt); r.UseDevice(true);
708  d.SetSize(width, mt); d.UseDevice(true);
709  z.SetSize(width, mt); z.UseDevice(true);
710 }
711 
712 void CGSolver::Mult(const Vector &b, Vector &x) const
713 {
714  int i;
715  double r0, den, nom, nom0, betanom, alpha, beta;
716 
717  x.UseDevice(true);
718  if (iterative_mode)
719  {
720  oper->Mult(x, r);
721  subtract(b, r, r); // r = b - A x
722  }
723  else
724  {
725  r = b;
726  x = 0.0;
727  }
728 
729  if (prec)
730  {
731  prec->Mult(r, z); // z = B r
732  d = z;
733  }
734  else
735  {
736  d = r;
737  }
738  nom0 = nom = Dot(d, r);
739  if (nom0 >= 0.0) { initial_norm = sqrt(nom0); }
740  MFEM_ASSERT(IsFinite(nom), "nom = " << nom);
742  {
743  mfem::out << " Iteration : " << setw(3) << 0 << " (B r, r) = "
744  << nom << (print_options.first_and_last ? " ...\n" : "\n");
745  }
746  Monitor(0, nom, r, x);
747 
748  if (nom < 0.0)
749  {
751  {
752  mfem::out << "PCG: The preconditioner is not positive definite. (Br, r) = "
753  << nom << '\n';
754  }
755  converged = false;
756  final_iter = 0;
757  initial_norm = nom;
758  final_norm = nom;
759  return;
760  }
761  r0 = std::max(nom*rel_tol*rel_tol, abs_tol*abs_tol);
762  if (nom <= r0)
763  {
764  converged = true;
765  final_iter = 0;
766  final_norm = sqrt(nom);
767  return;
768  }
769 
770  oper->Mult(d, z); // z = A d
771  den = Dot(z, d);
772  MFEM_ASSERT(IsFinite(den), "den = " << den);
773  if (den <= 0.0)
774  {
775  if (Dot(d, d) > 0.0 && print_options.warnings)
776  {
777  mfem::out << "PCG: The operator is not positive definite. (Ad, d) = "
778  << den << '\n';
779  }
780  if (den == 0.0)
781  {
782  converged = false;
783  final_iter = 0;
784  final_norm = sqrt(nom);
785  return;
786  }
787  }
788 
789  // start iteration
790  converged = false;
792  for (i = 1; true; )
793  {
794  alpha = nom/den;
795  add(x, alpha, d, x); // x = x + alpha d
796  add(r, -alpha, z, r); // r = r - alpha A d
797 
798  if (prec)
799  {
800  prec->Mult(r, z); // z = B r
801  betanom = Dot(r, z);
802  }
803  else
804  {
805  betanom = Dot(r, r);
806  }
807  MFEM_ASSERT(IsFinite(betanom), "betanom = " << betanom);
808  if (betanom < 0.0)
809  {
811  {
812  mfem::out << "PCG: The preconditioner is not positive definite. (Br, r) = "
813  << betanom << '\n';
814  }
815  converged = false;
816  final_iter = i;
817  break;
818  }
819 
821  {
822  mfem::out << " Iteration : " << setw(3) << i << " (B r, r) = "
823  << betanom << std::endl;
824  }
825 
826  Monitor(i, betanom, r, x);
827 
828  if (betanom <= r0)
829  {
830  converged = true;
831  final_iter = i;
832  break;
833  }
834 
835  if (++i > max_iter)
836  {
837  break;
838  }
839 
840  beta = betanom/nom;
841  if (prec)
842  {
843  add(z, beta, d, d); // d = z + beta d
844  }
845  else
846  {
847  add(r, beta, d, d);
848  }
849  oper->Mult(d, z); // z = A d
850  den = Dot(d, z);
851  MFEM_ASSERT(IsFinite(den), "den = " << den);
852  if (den <= 0.0)
853  {
854  if (Dot(d, d) > 0.0 && print_options.warnings)
855  {
856  mfem::out << "PCG: The operator is not positive definite. (Ad, d) = "
857  << den << '\n';
858  }
859  if (den == 0.0)
860  {
861  final_iter = i;
862  break;
863  }
864  }
865  nom = betanom;
866  }
868  {
869  mfem::out << " Iteration : " << setw(3) << final_iter << " (B r, r) = "
870  << betanom << '\n';
871  }
873  {
874  mfem::out << "PCG: Number of iterations: " << final_iter << '\n';
875  }
878  {
879  const auto arf = pow (betanom/nom0, 0.5/final_iter);
880  mfem::out << "Average reduction factor = " << arf << '\n';
881  }
883  {
884  mfem::out << "PCG: No convergence!" << '\n';
885  }
886 
887  final_norm = sqrt(betanom);
888 
889  Monitor(final_iter, final_norm, r, x, true);
890 }
891 
892 void CG(const Operator &A, const Vector &b, Vector &x,
893  int print_iter, int max_num_iter,
894  double RTOLERANCE, double ATOLERANCE)
895 {
896  MFEM_PERF_FUNCTION;
897 
898  CGSolver cg;
899  cg.SetPrintLevel(print_iter);
900  cg.SetMaxIter(max_num_iter);
901  cg.SetRelTol(sqrt(RTOLERANCE));
902  cg.SetAbsTol(sqrt(ATOLERANCE));
903  cg.SetOperator(A);
904  cg.Mult(b, x);
905 }
906 
907 void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x,
908  int print_iter, int max_num_iter,
909  double RTOLERANCE, double ATOLERANCE)
910 {
911  MFEM_PERF_FUNCTION;
912 
913  CGSolver pcg;
914  pcg.SetPrintLevel(print_iter);
915  pcg.SetMaxIter(max_num_iter);
916  pcg.SetRelTol(sqrt(RTOLERANCE));
917  pcg.SetAbsTol(sqrt(ATOLERANCE));
918  pcg.SetOperator(A);
919  pcg.SetPreconditioner(B);
920  pcg.Mult(b, x);
921 }
922 
923 
924 inline void GeneratePlaneRotation(double &dx, double &dy,
925  double &cs, double &sn)
926 {
927  if (dy == 0.0)
928  {
929  cs = 1.0;
930  sn = 0.0;
931  }
932  else if (fabs(dy) > fabs(dx))
933  {
934  double temp = dx / dy;
935  sn = 1.0 / sqrt( 1.0 + temp*temp );
936  cs = temp * sn;
937  }
938  else
939  {
940  double temp = dy / dx;
941  cs = 1.0 / sqrt( 1.0 + temp*temp );
942  sn = temp * cs;
943  }
944 }
945 
946 inline void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
947 {
948  double temp = cs * dx + sn * dy;
949  dy = -sn * dx + cs * dy;
950  dx = temp;
951 }
952 
953 inline void Update(Vector &x, int k, DenseMatrix &h, Vector &s,
954  Array<Vector*> &v)
955 {
956  Vector y(s);
957 
958  // Backsolve:
959  for (int i = k; i >= 0; i--)
960  {
961  y(i) /= h(i,i);
962  for (int j = i - 1; j >= 0; j--)
963  {
964  y(j) -= h(j,i) * y(i);
965  }
966  }
967 
968  for (int j = 0; j <= k; j++)
969  {
970  x.Add(y(j), *v[j]);
971  }
972 }
973 
974 void GMRESSolver::Mult(const Vector &b, Vector &x) const
975 {
976  // Generalized Minimum Residual method following the algorithm
977  // on p. 20 of the SIAM Templates book.
978 
979  int n = width;
980 
981  DenseMatrix H(m+1, m);
982  Vector s(m+1), cs(m+1), sn(m+1);
983  Vector r(n), w(n);
984  Array<Vector *> v;
985 
986  double resid;
987  int i, j, k;
988 
989  if (iterative_mode)
990  {
991  oper->Mult(x, r);
992  }
993  else
994  {
995  x = 0.0;
996  }
997 
998  if (prec)
999  {
1000  if (iterative_mode)
1001  {
1002  subtract(b, r, w);
1003  prec->Mult(w, r); // r = M (b - A x)
1004  }
1005  else
1006  {
1007  prec->Mult(b, r);
1008  }
1009  }
1010  else
1011  {
1012  if (iterative_mode)
1013  {
1014  subtract(b, r, r);
1015  }
1016  else
1017  {
1018  r = b;
1019  }
1020  }
1021  double beta = initial_norm = Norm(r); // beta = ||r||
1022  MFEM_ASSERT(IsFinite(beta), "beta = " << beta);
1023 
1024  final_norm = std::max(rel_tol*beta, abs_tol);
1025 
1026  if (beta <= final_norm)
1027  {
1028  final_norm = beta;
1029  final_iter = 0;
1030  converged = true;
1031  goto finish;
1032  }
1033 
1035  {
1036  mfem::out << " Pass : " << setw(2) << 1
1037  << " Iteration : " << setw(3) << 0
1038  << " ||B r|| = " << beta
1039  << (print_options.first_and_last ? " ...\n" : "\n");
1040  }
1041 
1042  Monitor(0, beta, r, x);
1043 
1044  v.SetSize(m+1, NULL);
1045 
1046  for (j = 1; j <= max_iter; )
1047  {
1048  if (v[0] == NULL) { v[0] = new Vector(n); }
1049  v[0]->Set(1.0/beta, r);
1050  s = 0.0; s(0) = beta;
1051 
1052  for (i = 0; i < m && j <= max_iter; i++, j++)
1053  {
1054  if (prec)
1055  {
1056  oper->Mult(*v[i], r);
1057  prec->Mult(r, w); // w = M A v[i]
1058  }
1059  else
1060  {
1061  oper->Mult(*v[i], w);
1062  }
1063 
1064  for (k = 0; k <= i; k++)
1065  {
1066  H(k,i) = Dot(w, *v[k]); // H(k,i) = w * v[k]
1067  w.Add(-H(k,i), *v[k]); // w -= H(k,i) * v[k]
1068  }
1069 
1070  H(i+1,i) = Norm(w); // H(i+1,i) = ||w||
1071  MFEM_ASSERT(IsFinite(H(i+1,i)), "Norm(w) = " << H(i+1,i));
1072  if (v[i+1] == NULL) { v[i+1] = new Vector(n); }
1073  v[i+1]->Set(1.0/H(i+1,i), w); // v[i+1] = w / H(i+1,i)
1074 
1075  for (k = 0; k < i; k++)
1076  {
1077  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
1078  }
1079 
1080  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1081  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1082  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
1083 
1084  resid = fabs(s(i+1));
1085  MFEM_ASSERT(IsFinite(resid), "resid = " << resid);
1086 
1087  if (resid <= final_norm)
1088  {
1089  Update(x, i, H, s, v);
1090  final_norm = resid;
1091  final_iter = j;
1092  converged = true;
1093  goto finish;
1094  }
1095 
1097  {
1098  mfem::out << " Pass : " << setw(2) << (j-1)/m+1
1099  << " Iteration : " << setw(3) << j
1100  << " ||B r|| = " << resid << '\n';
1101  }
1102 
1103  Monitor(j, resid, r, x);
1104  }
1105 
1106  if (print_options.iterations && j <= max_iter)
1107  {
1108  mfem::out << "Restarting..." << '\n';
1109  }
1110 
1111  Update(x, i-1, H, s, v);
1112 
1113  oper->Mult(x, r);
1114  if (prec)
1115  {
1116  subtract(b, r, w);
1117  prec->Mult(w, r); // r = M (b - A x)
1118  }
1119  else
1120  {
1121  subtract(b, r, r);
1122  }
1123  beta = Norm(r); // beta = ||r||
1124  MFEM_ASSERT(IsFinite(beta), "beta = " << beta);
1125  if (beta <= final_norm)
1126  {
1127  final_norm = beta;
1128  final_iter = j;
1129  converged = true;
1130  goto finish;
1131  }
1132  }
1133 
1134  final_norm = beta;
1135  final_iter = max_iter;
1136  converged = false;
1137 
1138 finish:
1140  {
1141  mfem::out << " Pass : " << setw(2) << (j-1)/m+1
1142  << " Iteration : " << setw(3) << final_iter
1143  << " ||B r|| = " << resid << '\n';
1144  }
1146  {
1147  mfem::out << "GMRES: Number of iterations: " << final_iter << '\n';
1148  }
1150  {
1151  mfem::out << "GMRES: No convergence!\n";
1152  }
1153 
1154  Monitor(final_iter, final_norm, r, x, true);
1155 
1156  for (i = 0; i < v.Size(); i++)
1157  {
1158  delete v[i];
1159  }
1160 }
1161 
1162 void FGMRESSolver::Mult(const Vector &b, Vector &x) const
1163 {
1164  DenseMatrix H(m+1,m);
1165  Vector s(m+1), cs(m+1), sn(m+1);
1166  Vector r(b.Size());
1167 
1168  int i, j, k;
1169 
1170 
1171  if (iterative_mode)
1172  {
1173  oper->Mult(x, r);
1174  subtract(b,r,r);
1175  }
1176  else
1177  {
1178  x = 0.;
1179  r = b;
1180  }
1181  double beta = initial_norm = Norm(r); // beta = ||r||
1182  // We need to preallocate this to report the correct result in the case of
1183  // no convergence.
1184  double resid;
1185  MFEM_ASSERT(IsFinite(beta), "beta = " << beta);
1186 
1187  final_norm = std::max(rel_tol*beta, abs_tol);
1188 
1189  if (beta <= final_norm)
1190  {
1191  final_norm = beta;
1192  final_iter = 0;
1193  converged = true;
1194  return;
1195  }
1196 
1198  {
1199  mfem::out << " Pass : " << setw(2) << 1
1200  << " Iteration : " << setw(3) << 0
1201  << " || r || = " << beta
1202  << (print_options.first_and_last ? " ...\n" : "\n");
1203  }
1204 
1205  Monitor(0, beta, r, x);
1206 
1207  Array<Vector*> v(m+1);
1208  Array<Vector*> z(m+1);
1209  for (i= 0; i<=m; i++)
1210  {
1211  v[i] = NULL;
1212  z[i] = NULL;
1213  }
1214 
1215  j = 1;
1216  while (j <= max_iter)
1217  {
1218  if (v[0] == NULL) { v[0] = new Vector(b.Size()); }
1219  (*v[0]) = 0.0;
1220  v[0] -> Add (1.0/beta, r); // v[0] = r / ||r||
1221  s = 0.0; s(0) = beta;
1222 
1223  for (i = 0; i < m && j <= max_iter; i++, j++)
1224  {
1225 
1226  if (z[i] == NULL) { z[i] = new Vector(b.Size()); }
1227  (*z[i]) = 0.0;
1228 
1229  if (prec)
1230  {
1231  prec->Mult(*v[i], *z[i]);
1232  }
1233  else
1234  {
1235  (*z[i]) = (*v[i]);
1236  }
1237  oper->Mult(*z[i], r);
1238 
1239  for (k = 0; k <= i; k++)
1240  {
1241  H(k,i) = Dot( r, *v[k]); // H(k,i) = r * v[k]
1242  r.Add(-H(k,i), (*v[k])); // r -= H(k,i) * v[k]
1243  }
1244 
1245  H(i+1,i) = Norm(r); // H(i+1,i) = ||r||
1246  if (v[i+1] == NULL) { v[i+1] = new Vector(b.Size()); }
1247  (*v[i+1]) = 0.0;
1248  v[i+1] -> Add (1.0/H(i+1,i), r); // v[i+1] = r / H(i+1,i)
1249 
1250  for (k = 0; k < i; k++)
1251  {
1252  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
1253  }
1254 
1255  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1256  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
1257  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
1258 
1259  resid = fabs(s(i+1));
1260  MFEM_ASSERT(IsFinite(resid), "resid = " << resid);
1262  resid <= final_norm))
1263  {
1264  mfem::out << " Pass : " << setw(2) << (j-1)/m+1
1265  << " Iteration : " << setw(3) << j
1266  << " || r || = " << resid << endl;
1267  }
1268  Monitor(j, resid, r, x, resid <= final_norm);
1269 
1270  if (resid <= final_norm)
1271  {
1272  Update(x, i, H, s, z);
1273  final_norm = resid;
1274  final_iter = j;
1275  converged = true;
1276 
1277  if (print_options.summary)
1278  {
1279  mfem::out << "FGMRES: Number of iterations: " << final_iter << '\n';
1280  }
1281 
1282  for (i= 0; i<=m; i++)
1283  {
1284  if (v[i]) { delete v[i]; }
1285  if (z[i]) { delete z[i]; }
1286  }
1287  return;
1288  }
1289  }
1290 
1292  {
1293  mfem::out << "Restarting..." << endl;
1294  }
1295 
1296  Update(x, i-1, H, s, z);
1297 
1298  oper->Mult(x, r);
1299  subtract(b,r,r);
1300  beta = Norm(r);
1301  MFEM_ASSERT(IsFinite(beta), "beta = " << beta);
1302  if (beta <= final_norm)
1303  {
1304  final_norm = beta;
1305  final_iter = j;
1306  converged = true;
1307 
1308  break;
1309  }
1310  }
1311 
1312  // Clean buffers up
1313  for (i = 0; i <= m; i++)
1314  {
1315  if (v[i]) { delete v[i]; }
1316  if (z[i]) { delete z[i]; }
1317  }
1318  converged = false;
1319 
1320  // Note: j is off by one when we arrive here
1322  {
1323  mfem::out << " Pass : " << setw(2) << (j-1)/m+1
1324  << " Iteration : " << setw(3) << j-1
1325  << " || r || = " << resid << endl;
1326  }
1328  {
1329  mfem::out << "FGMRES: Number of iterations: " << j-1 << '\n';
1330  }
1332  {
1333  mfem::out << "FGMRES: No convergence!\n";
1334  }
1335 }
1336 
1337 
1338 int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M,
1339  int &max_iter, int m, double &tol, double atol, int printit)
1340 {
1341  MFEM_PERF_FUNCTION;
1342 
1343  GMRESSolver gmres;
1344  gmres.SetPrintLevel(printit);
1345  gmres.SetMaxIter(max_iter);
1346  gmres.SetKDim(m);
1347  gmres.SetRelTol(sqrt(tol));
1348  gmres.SetAbsTol(sqrt(atol));
1349  gmres.SetOperator(A);
1350  gmres.SetPreconditioner(M);
1351  gmres.Mult(b, x);
1352  max_iter = gmres.GetNumIterations();
1353  tol = gmres.GetFinalNorm()*gmres.GetFinalNorm();
1354  return gmres.GetConverged();
1355 }
1356 
1357 void GMRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
1358  int print_iter, int max_num_iter, int m, double rtol, double atol)
1359 {
1360  GMRES(A, x, b, B, max_num_iter, m, rtol, atol, print_iter);
1361 }
1362 
1363 
1365 {
1366  p.SetSize(width);
1367  phat.SetSize(width);
1368  s.SetSize(width);
1369  shat.SetSize(width);
1370  t.SetSize(width);
1371  v.SetSize(width);
1372  r.SetSize(width);
1373  rtilde.SetSize(width);
1374 }
1375 
1376 void BiCGSTABSolver::Mult(const Vector &b, Vector &x) const
1377 {
1378  // BiConjugate Gradient Stabilized method following the algorithm
1379  // on p. 27 of the SIAM Templates book.
1380 
1381  int i;
1382  double resid, tol_goal;
1383  double rho_1, rho_2=1.0, alpha=1.0, beta, omega=1.0;
1384 
1385  if (iterative_mode)
1386  {
1387  oper->Mult(x, r);
1388  subtract(b, r, r); // r = b - A x
1389  }
1390  else
1391  {
1392  x = 0.0;
1393  r = b;
1394  }
1395  rtilde = r;
1396 
1397  resid = initial_norm = Norm(r);
1398  MFEM_ASSERT(IsFinite(resid), "resid = " << resid);
1400  {
1401  mfem::out << " Iteration : " << setw(3) << 0
1402  << " ||r|| = " << resid << (print_options.first_and_last ? " ...\n" : "\n");
1403  }
1404 
1405  Monitor(0, resid, r, x);
1406 
1407  tol_goal = std::max(resid*rel_tol, abs_tol);
1408 
1409  if (resid <= tol_goal)
1410  {
1411  final_norm = resid;
1412  final_iter = 0;
1413  converged = true;
1414  return;
1415  }
1416 
1417  for (i = 1; i <= max_iter; i++)
1418  {
1419  rho_1 = Dot(rtilde, r);
1420  if (rho_1 == 0)
1421  {
1423  {
1424  mfem::out << " Iteration : " << setw(3) << i
1425  << " ||r|| = " << resid << '\n';
1426  }
1427 
1428  Monitor(i, resid, r, x);
1429 
1430  final_norm = resid;
1431  final_iter = i;
1432  converged = false;
1434  {
1435  mfem::out << "BiCGStab: Number of iterations: " << final_iter << '\n';
1436  }
1437  if (print_options.warnings)
1438  {
1439  mfem::out << "BiCGStab: No convergence!\n";
1440  }
1441  return;
1442  }
1443  if (i == 1)
1444  {
1445  p = r;
1446  }
1447  else
1448  {
1449  beta = (rho_1/rho_2) * (alpha/omega);
1450  add(p, -omega, v, p); // p = p - omega * v
1451  add(r, beta, p, p); // p = r + beta * p
1452  }
1453  if (prec)
1454  {
1455  prec->Mult(p, phat); // phat = M^{-1} * p
1456  }
1457  else
1458  {
1459  phat = p;
1460  }
1461  oper->Mult(phat, v); // v = A * phat
1462  alpha = rho_1 / Dot(rtilde, v);
1463  add(r, -alpha, v, s); // s = r - alpha * v
1464  resid = Norm(s);
1465  MFEM_ASSERT(IsFinite(resid), "resid = " << resid);
1466  if (resid < tol_goal)
1467  {
1468  x.Add(alpha, phat); // x = x + alpha * phat
1470  {
1471  mfem::out << " Iteration : " << setw(3) << i
1472  << " ||s|| = " << resid << '\n';
1473  }
1474  final_norm = resid;
1475  final_iter = i;
1476  converged = true;
1478  {
1479  mfem::out << "BiCGStab: Number of iterations: " << final_iter << '\n';
1480  }
1481  return;
1482  }
1484  {
1485  mfem::out << " Iteration : " << setw(3) << i
1486  << " ||s|| = " << resid;
1487  }
1488  Monitor(i, resid, r, x);
1489  if (prec)
1490  {
1491  prec->Mult(s, shat); // shat = M^{-1} * s
1492  }
1493  else
1494  {
1495  shat = s;
1496  }
1497  oper->Mult(shat, t); // t = A * shat
1498  omega = Dot(t, s) / Dot(t, t);
1499  x.Add(alpha, phat); // x += alpha * phat
1500  x.Add(omega, shat); // x += omega * shat
1501  add(s, -omega, t, r); // r = s - omega * t
1502 
1503  rho_2 = rho_1;
1504  resid = Norm(r);
1505  MFEM_ASSERT(IsFinite(resid), "resid = " << resid);
1507  {
1508  mfem::out << " ||r|| = " << resid << '\n';
1509  }
1510  Monitor(i, resid, r, x);
1511  if (resid < tol_goal)
1512  {
1513  final_norm = resid;
1514  final_iter = i;
1515  converged = true;
1517  {
1518  mfem::out << " Iteration : " << setw(3) << i
1519  << " ||r|| = " << resid << '\n';
1520  }
1522  {
1523  mfem::out << "BiCGStab: Number of iterations: " << final_iter << '\n';
1524  }
1525  return;
1526  }
1527  if (omega == 0)
1528  {
1529  final_norm = resid;
1530  final_iter = i;
1531  converged = false;
1533  {
1534  mfem::out << " Iteration : " << setw(3) << i
1535  << " ||r|| = " << resid << '\n';
1536  }
1538  {
1539  mfem::out << "BiCGStab: Number of iterations: " << final_iter << '\n';
1540  }
1541  if (print_options.warnings)
1542  {
1543  mfem::out << "BiCGStab: No convergence!\n";
1544  }
1545  return;
1546  }
1547  }
1548 
1549  final_norm = resid;
1550  final_iter = max_iter;
1551  converged = false;
1552 
1554  {
1555  mfem::out << " Iteration : " << setw(3) << final_iter
1556  << " ||r|| = " << resid << '\n';
1557  }
1559  {
1560  mfem::out << "BiCGStab: Number of iterations: " << final_iter << '\n';
1561  }
1562  if (print_options.warnings)
1563  {
1564  mfem::out << "BiCGStab: No convergence!\n";
1565  }
1566 }
1567 
1568 int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M,
1569  int &max_iter, double &tol, double atol, int printit)
1570 {
1571  BiCGSTABSolver bicgstab;
1572  bicgstab.SetPrintLevel(printit);
1573  bicgstab.SetMaxIter(max_iter);
1574  bicgstab.SetRelTol(sqrt(tol));
1575  bicgstab.SetAbsTol(sqrt(atol));
1576  bicgstab.SetOperator(A);
1577  bicgstab.SetPreconditioner(M);
1578  bicgstab.Mult(b, x);
1579  max_iter = bicgstab.GetNumIterations();
1580  tol = bicgstab.GetFinalNorm()*bicgstab.GetFinalNorm();
1581  return bicgstab.GetConverged();
1582 }
1583 
1584 void BiCGSTAB(const Operator &A, Solver &B, const Vector &b, Vector &x,
1585  int print_iter, int max_num_iter, double rtol, double atol)
1586 {
1587  BiCGSTAB(A, x, b, B, max_num_iter, rtol, atol, print_iter);
1588 }
1589 
1590 
1592 {
1594  v0.SetSize(width);
1595  v1.SetSize(width);
1596  w0.SetSize(width);
1597  w1.SetSize(width);
1598  q.SetSize(width);
1599  if (prec)
1600  {
1601  u1.SetSize(width);
1602  }
1603 
1604  v0.UseDevice(true);
1605  v1.UseDevice(true);
1606  w0.UseDevice(true);
1607  w1.UseDevice(true);
1608  q.UseDevice(true);
1609  u1.UseDevice(true);
1610 }
1611 
1612 void MINRESSolver::Mult(const Vector &b, Vector &x) const
1613 {
1614  // Based on the MINRES algorithm on p. 86, Fig. 6.9 in
1615  // "Iterative Krylov Methods for Large Linear Systems",
1616  // by Henk A. van der Vorst, 2003.
1617  // Extended to support an SPD preconditioner.
1618 
1619  b.UseDevice(true);
1620  x.UseDevice(true);
1621 
1622  int it;
1623  double beta, eta, gamma0, gamma1, sigma0, sigma1;
1624  double alpha, delta, rho1, rho2, rho3, norm_goal;
1625  Vector *z = (prec) ? &u1 : &v1;
1626 
1627  converged = true;
1628 
1629  if (!iterative_mode)
1630  {
1631  v1 = b;
1632  x = 0.;
1633  }
1634  else
1635  {
1636  oper->Mult(x, v1);
1637  subtract(b, v1, v1);
1638  }
1639 
1640  if (prec)
1641  {
1642  prec->Mult(v1, u1);
1643  }
1644  eta = beta = initial_norm = sqrt(Dot(*z, v1));
1645  MFEM_ASSERT(IsFinite(eta), "eta = " << eta);
1646  gamma0 = gamma1 = 1.;
1647  sigma0 = sigma1 = 0.;
1648 
1649  norm_goal = std::max(rel_tol*eta, abs_tol);
1650 
1651  if (eta <= norm_goal)
1652  {
1653  it = 0;
1654  goto loop_end;
1655  }
1656 
1658  {
1659  mfem::out << "MINRES: iteration " << setw(3) << 0 << ": ||r||_B = "
1660  << eta << (print_options.first_and_last ? " ..." : "") << '\n';
1661  }
1662  Monitor(0, eta, *z, x);
1663 
1664  for (it = 1; it <= max_iter; it++)
1665  {
1666  v1 /= beta;
1667  if (prec)
1668  {
1669  u1 /= beta;
1670  }
1671  oper->Mult(*z, q);
1672  alpha = Dot(*z, q);
1673  MFEM_ASSERT(IsFinite(alpha), "alpha = " << alpha);
1674  if (it > 1) // (v0 == 0) for (it == 1)
1675  {
1676  q.Add(-beta, v0);
1677  }
1678  add(q, -alpha, v1, v0);
1679 
1680  delta = gamma1*alpha - gamma0*sigma1*beta;
1681  rho3 = sigma0*beta;
1682  rho2 = sigma1*alpha + gamma0*gamma1*beta;
1683  if (!prec)
1684  {
1685  beta = Norm(v0);
1686  }
1687  else
1688  {
1689  prec->Mult(v0, q);
1690  beta = sqrt(Dot(v0, q));
1691  }
1692  MFEM_ASSERT(IsFinite(beta), "beta = " << beta);
1693  rho1 = hypot(delta, beta);
1694 
1695  if (it == 1)
1696  {
1697  w0.Set(1./rho1, *z); // (w0 == 0) and (w1 == 0)
1698  }
1699  else if (it == 2)
1700  {
1701  add(1./rho1, *z, -rho2/rho1, w1, w0); // (w0 == 0)
1702  }
1703  else
1704  {
1705  add(-rho3/rho1, w0, -rho2/rho1, w1, w0);
1706  w0.Add(1./rho1, *z);
1707  }
1708 
1709  gamma0 = gamma1;
1710  gamma1 = delta/rho1;
1711 
1712  x.Add(gamma1*eta, w0);
1713 
1714  sigma0 = sigma1;
1715  sigma1 = beta/rho1;
1716 
1717  eta = -sigma1*eta;
1718  MFEM_ASSERT(IsFinite(eta), "eta = " << eta);
1719 
1720  if (fabs(eta) <= norm_goal)
1721  {
1722  goto loop_end;
1723  }
1724 
1726  {
1727  mfem::out << "MINRES: iteration " << setw(3) << it << ": ||r||_B = "
1728  << fabs(eta) << '\n';
1729  }
1730  Monitor(it, fabs(eta), *z, x);
1731 
1732  if (prec)
1733  {
1734  Swap(u1, q);
1735  }
1736  Swap(v0, v1);
1737  Swap(w0, w1);
1738  }
1739  converged = false;
1740  it--;
1741 
1742 loop_end:
1743  final_iter = it;
1744  final_norm = fabs(eta);
1745 
1747  {
1748  mfem::out << "MINRES: iteration " << setw(3) << it << ": ||r||_B = "
1749  << fabs(eta) << '\n';
1750  }
1751 
1753  {
1754  mfem::out << "MINRES: Number of iterations: " << setw(3) << final_iter << '\n';
1755  }
1756 
1757  Monitor(final_iter, final_norm, *z, x, true);
1758 
1759  // if (print_options.iteration_details || (!converged && print_options.errors))
1760  // {
1761  // oper->Mult(x, v1);
1762  // subtract(b, v1, v1);
1763  // if (prec)
1764  // {
1765  // prec->Mult(v1, u1);
1766  // }
1767  // eta = sqrt(Dot(*z, v1));
1768  // mfem::out << "MINRES: iteration " << setw(3) << it << '\n'
1769  // << " ||r||_B = " << eta << " (re-computed)" << '\n';
1770  // }
1771 
1772  if (!converged && (print_options.warnings))
1773  {
1774  mfem::out << "MINRES: No convergence!\n";
1775  }
1776 }
1777 
1778 void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it,
1779  int max_it, double rtol, double atol)
1780 {
1781  MFEM_PERF_FUNCTION;
1782 
1783  MINRESSolver minres;
1784  minres.SetPrintLevel(print_it);
1785  minres.SetMaxIter(max_it);
1786  minres.SetRelTol(sqrt(rtol));
1787  minres.SetAbsTol(sqrt(atol));
1788  minres.SetOperator(A);
1789  minres.Mult(b, x);
1790 }
1791 
1792 void MINRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
1793  int print_it, int max_it, double rtol, double atol)
1794 {
1795  MINRESSolver minres;
1796  minres.SetPrintLevel(print_it);
1797  minres.SetMaxIter(max_it);
1798  minres.SetRelTol(sqrt(rtol));
1799  minres.SetAbsTol(sqrt(atol));
1800  minres.SetOperator(A);
1801  minres.SetPreconditioner(B);
1802  minres.Mult(b, x);
1803 }
1804 
1805 
1807 {
1808  oper = &op;
1809  height = op.Height();
1810  width = op.Width();
1811  MFEM_ASSERT(height == width, "square Operator is required.");
1812 
1813  r.SetSize(width);
1814  c.SetSize(width);
1815 }
1816 
1817 void NewtonSolver::Mult(const Vector &b, Vector &x) const
1818 {
1819  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
1820  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
1821 
1822  int it;
1823  double norm0, norm, norm_goal;
1824  const bool have_b = (b.Size() == Height());
1825 
1826  if (!iterative_mode)
1827  {
1828  x = 0.0;
1829  }
1830 
1831  ProcessNewState(x);
1832 
1833  oper->Mult(x, r);
1834  if (have_b)
1835  {
1836  r -= b;
1837  }
1838 
1839  norm0 = norm = initial_norm = Norm(r);
1841  {
1842  mfem::out << "Newton iteration " << setw(2) << 0
1843  << " : ||r|| = " << norm << "...\n";
1844  }
1845  norm_goal = std::max(rel_tol*norm, abs_tol);
1846 
1847  prec->iterative_mode = false;
1848 
1849  // x_{i+1} = x_i - [DF(x_i)]^{-1} [F(x_i)-b]
1850  for (it = 0; true; it++)
1851  {
1852  MFEM_ASSERT(IsFinite(norm), "norm = " << norm);
1854  {
1855  mfem::out << "Newton iteration " << setw(2) << it
1856  << " : ||r|| = " << norm;
1857  if (it > 0)
1858  {
1859  mfem::out << ", ||r||/||r_0|| = " << norm/norm0;
1860  }
1861  mfem::out << '\n';
1862  }
1863  Monitor(it, norm, r, x);
1864 
1865  if (norm <= norm_goal)
1866  {
1867  converged = true;
1868  break;
1869  }
1870 
1871  if (it >= max_iter)
1872  {
1873  converged = false;
1874  break;
1875  }
1876 
1877  grad = &oper->GetGradient(x);
1878  prec->SetOperator(*grad);
1879 
1880  if (lin_rtol_type)
1881  {
1882  AdaptiveLinRtolPreSolve(x, it, norm);
1883  }
1884 
1885  prec->Mult(r, c); // c = [DF(x_i)]^{-1} [F(x_i)-b]
1886 
1887  if (lin_rtol_type)
1888  {
1889  AdaptiveLinRtolPostSolve(c, r, it, norm);
1890  }
1891 
1892  const double c_scale = ComputeScalingFactor(x, b);
1893  if (c_scale == 0.0)
1894  {
1895  converged = false;
1896  break;
1897  }
1898  add(x, -c_scale, c, x);
1899 
1900  ProcessNewState(x);
1901 
1902  oper->Mult(x, r);
1903  if (have_b)
1904  {
1905  r -= b;
1906  }
1907  norm = Norm(r);
1908  }
1909 
1910  final_iter = it;
1911  final_norm = norm;
1912 
1915  {
1916  mfem::out << "Newton: Number of iterations: " << final_iter << '\n'
1917  << " ||r|| = " << final_norm << '\n';
1918  }
1920  {
1921  mfem::out << "Newton: No convergence!\n";
1922  }
1923 }
1924 
1926  const double rtol0,
1927  const double rtol_max,
1928  const double alpha_,
1929  const double gamma_)
1930 {
1931  lin_rtol_type = type;
1932  lin_rtol0 = rtol0;
1933  lin_rtol_max = rtol_max;
1934  this->alpha = alpha_;
1935  this->gamma = gamma_;
1936 }
1937 
1939  const int it,
1940  const double fnorm) const
1941 {
1942  // Assume that when adaptive linear solver relative tolerance is activated,
1943  // we are working with an iterative solver.
1944  auto iterative_solver = static_cast<IterativeSolver *>(prec);
1945  // Adaptive linear solver relative tolerance
1946  double eta;
1947  // Safeguard threshold
1948  double sg_threshold = 0.1;
1949 
1950  if (it == 0)
1951  {
1952  eta = lin_rtol0;
1953  }
1954  else
1955  {
1956  if (lin_rtol_type == 1)
1957  {
1958  // eta = gamma * abs(||F(x1)|| - ||F(x0) + DF(x0) s0||) / ||F(x0)||
1959  eta = gamma * abs(fnorm - lnorm_last) / fnorm_last;
1960  }
1961  else if (lin_rtol_type == 2)
1962  {
1963  // eta = gamma * (||F(x1)|| / ||F(x0)||)^alpha
1964  eta = gamma * pow(fnorm / fnorm_last, alpha);
1965  }
1966  else
1967  {
1968  MFEM_ABORT("Unknown adaptive linear solver rtol version");
1969  }
1970 
1971  // Safeguard rtol from "oversolving" ?!
1972  const double sg_eta = gamma * pow(eta_last, alpha);
1973  if (sg_eta > sg_threshold) { eta = std::max(eta, sg_eta); }
1974  }
1975 
1976  eta = std::min(eta, lin_rtol_max);
1977  iterative_solver->SetRelTol(eta);
1978  eta_last = eta;
1980  {
1981  mfem::out << "Eisenstat-Walker rtol = " << eta << "\n";
1982  }
1983 }
1984 
1986  const Vector &b,
1987  const int it,
1988  const double fnorm) const
1989 {
1990  fnorm_last = fnorm;
1991 
1992  // If version 1 is chosen, the true linear residual norm has to be computed
1993  // and in most cases we can only retrieve the preconditioned linear residual
1994  // norm.
1995  if (lin_rtol_type == 1)
1996  {
1997  // lnorm_last = ||F(x0) + DF(x0) s0||
1998  Vector linres(x.Size());
1999  grad->Mult(x, linres);
2000  linres -= b;
2001  lnorm_last = Norm(linres);
2002  }
2003 }
2004 
2005 void LBFGSSolver::Mult(const Vector &b, Vector &x) const
2006 {
2007  MFEM_VERIFY(oper != NULL, "the Operator is not set (use SetOperator).");
2008 
2009  // Quadrature points that are checked for negative Jacobians etc.
2010  Vector sk, rk, yk, rho, alpha;
2011 
2012  // r - r_{k+1}, c - descent direction
2013  sk.SetSize(width); // x_{k+1}-x_k
2014  rk.SetSize(width); // nabla(f(x_{k}))
2015  yk.SetSize(width); // r_{k+1}-r_{k}
2016  rho.SetSize(m); // 1/(dot(yk,sk)
2017  alpha.SetSize(m); // rhok*sk'*c
2018  int last_saved_id = -1;
2019 
2020  int it;
2021  double norm0, norm, norm_goal;
2022  const bool have_b = (b.Size() == Height());
2023 
2024  if (!iterative_mode)
2025  {
2026  x = 0.0;
2027  }
2028 
2029  ProcessNewState(x);
2030 
2031  // r = F(x)-b
2032  oper->Mult(x, r);
2033  if (have_b) { r -= b; }
2034 
2035  c = r; // initial descent direction
2036 
2037  norm0 = norm = initial_norm = Norm(r);
2039  {
2040  mfem::out << "LBFGS iteration " << setw(2) << 0
2041  << " : ||r|| = " << norm << "...\n";
2042  }
2043  norm_goal = std::max(rel_tol*norm, abs_tol);
2044  for (it = 0; true; it++)
2045  {
2046  MFEM_ASSERT(IsFinite(norm), "norm = " << norm);
2048  {
2049  mfem::out << "LBFGS iteration " << it
2050  << " : ||r|| = " << norm;
2051  if (it > 0)
2052  {
2053  mfem::out << ", ||r||/||r_0|| = " << norm/norm0;
2054  }
2055  mfem::out << '\n';
2056  }
2057 
2058  if (norm <= norm_goal)
2059  {
2060  converged = true;
2061  break;
2062  }
2063 
2064  if (it >= max_iter)
2065  {
2066  converged = false;
2067  break;
2068  }
2069 
2070  rk = r;
2071  const double c_scale = ComputeScalingFactor(x, b);
2072  if (c_scale == 0.0)
2073  {
2074  converged = false;
2075  break;
2076  }
2077  add(x, -c_scale, c, x); // x_{k+1} = x_k - c_scale*c
2078 
2079  ProcessNewState(x);
2080 
2081  oper->Mult(x, r);
2082  if (have_b)
2083  {
2084  r -= b;
2085  }
2086 
2087  // LBFGS - construct descent direction
2088  subtract(r, rk, yk); // yk = r_{k+1} - r_{k}
2089  sk = c; sk *= -c_scale; //sk = x_{k+1} - x_{k} = -c_scale*c
2090  const double gamma = Dot(sk, yk)/Dot(yk, yk);
2091 
2092  // Save last m vectors
2093  last_saved_id = (last_saved_id == m-1) ? 0 : last_saved_id+1;
2094  *skArray[last_saved_id] = sk;
2095  *ykArray[last_saved_id] = yk;
2096 
2097  c = r;
2098  for (int i = last_saved_id; i > -1; i--)
2099  {
2100  rho(i) = 1.0/Dot((*skArray[i]),(*ykArray[i]));
2101  alpha(i) = rho(i)*Dot((*skArray[i]),c);
2102  add(c, -alpha(i), (*ykArray[i]), c);
2103  }
2104  if (it > m-1)
2105  {
2106  for (int i = m-1; i > last_saved_id; i--)
2107  {
2108  rho(i) = 1./Dot((*skArray[i]), (*ykArray[i]));
2109  alpha(i) = rho(i)*Dot((*skArray[i]),c);
2110  add(c, -alpha(i), (*ykArray[i]), c);
2111  }
2112  }
2113 
2114  c *= gamma; // scale search direction
2115  if (it > m-1)
2116  {
2117  for (int i = last_saved_id+1; i < m ; i++)
2118  {
2119  double betai = rho(i)*Dot((*ykArray[i]), c);
2120  add(c, alpha(i)-betai, (*skArray[i]), c);
2121  }
2122  }
2123  for (int i = 0; i < last_saved_id+1 ; i++)
2124  {
2125  double betai = rho(i)*Dot((*ykArray[i]), c);
2126  add(c, alpha(i)-betai, (*skArray[i]), c);
2127  }
2128 
2129  norm = Norm(r);
2130  }
2131 
2132  final_iter = it;
2133  final_norm = norm;
2134 
2137  {
2138  mfem::out << "LBFGS: Number of iterations: " << final_iter << '\n'
2139  << " ||r|| = " << final_norm << '\n';
2140  }
2142  {
2143  mfem::out << "LBFGS: No convergence!\n";
2144  }
2145 }
2146 
2147 int aGMRES(const Operator &A, Vector &x, const Vector &b,
2148  const Operator &M, int &max_iter,
2149  int m_max, int m_min, int m_step, double cf,
2150  double &tol, double &atol, int printit)
2151 {
2152  int n = A.Width();
2153 
2154  int m = m_max;
2155 
2156  DenseMatrix H(m+1,m);
2157  Vector s(m+1), cs(m+1), sn(m+1);
2158  Vector w(n), av(n);
2159 
2160  double r1, resid;
2161  int i, j, k;
2162 
2163  M.Mult(b,w);
2164  double normb = w.Norml2(); // normb = ||M b||
2165  if (normb == 0.0)
2166  {
2167  normb = 1;
2168  }
2169 
2170  Vector r(n);
2171  A.Mult(x, r);
2172  subtract(b,r,w);
2173  M.Mult(w, r); // r = M (b - A x)
2174  double beta = r.Norml2(); // beta = ||r||
2175 
2176  resid = beta / normb;
2177 
2178  if (resid * resid <= tol)
2179  {
2180  tol = resid * resid;
2181  max_iter = 0;
2182  return 0;
2183  }
2184 
2185  if (printit)
2186  {
2187  mfem::out << " Pass : " << setw(2) << 1
2188  << " Iteration : " << setw(3) << 0
2189  << " (r, r) = " << beta*beta << '\n';
2190  }
2191 
2192  tol *= (normb*normb);
2193  tol = (atol > tol) ? atol : tol;
2194 
2195  m = m_max;
2196  Array<Vector *> v(m+1);
2197  for (i= 0; i<=m; i++)
2198  {
2199  v[i] = new Vector(n);
2200  (*v[i]) = 0.0;
2201  }
2202 
2203  j = 1;
2204  while (j <= max_iter)
2205  {
2206  (*v[0]) = 0.0;
2207  v[0] -> Add (1.0/beta, r); // v[0] = r / ||r||
2208  s = 0.0; s(0) = beta;
2209 
2210  r1 = beta;
2211 
2212  for (i = 0; i < m && j <= max_iter; i++)
2213  {
2214  A.Mult((*v[i]),av);
2215  M.Mult(av,w); // w = M A v[i]
2216 
2217  for (k = 0; k <= i; k++)
2218  {
2219  H(k,i) = w * (*v[k]); // H(k,i) = w * v[k]
2220  w.Add(-H(k,i), (*v[k])); // w -= H(k,i) * v[k]
2221  }
2222 
2223  H(i+1,i) = w.Norml2(); // H(i+1,i) = ||w||
2224  (*v[i+1]) = 0.0;
2225  v[i+1] -> Add (1.0/H(i+1,i), w); // v[i+1] = w / H(i+1,i)
2226 
2227  for (k = 0; k < i; k++)
2228  {
2229  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
2230  }
2231 
2232  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
2233  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
2234  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
2235 
2236  resid = fabs(s(i+1));
2237  if (printit)
2238  {
2239  mfem::out << " Pass : " << setw(2) << j
2240  << " Iteration : " << setw(3) << i+1
2241  << " (r, r) = " << resid*resid << '\n';
2242  }
2243 
2244  if ( resid*resid < tol)
2245  {
2246  Update(x, i, H, s, v);
2247  tol = resid * resid;
2248  max_iter = j;
2249  for (i= 0; i<=m; i++)
2250  {
2251  delete v[i];
2252  }
2253  return 0;
2254  }
2255  }
2256 
2257  if (printit)
2258  {
2259  mfem::out << "Restarting..." << '\n';
2260  }
2261 
2262  Update(x, i-1, H, s, v);
2263 
2264  A.Mult(x, r);
2265  subtract(b,r,w);
2266  M.Mult(w, r); // r = M (b - A x)
2267  beta = r.Norml2(); // beta = ||r||
2268  if ( resid*resid < tol)
2269  {
2270  tol = resid * resid;
2271  max_iter = j;
2272  for (i= 0; i<=m; i++)
2273  {
2274  delete v[i];
2275  }
2276  return 0;
2277  }
2278 
2279  if (beta/r1 > cf)
2280  {
2281  if (m - m_step >= m_min)
2282  {
2283  m -= m_step;
2284  }
2285  else
2286  {
2287  m = m_max;
2288  }
2289  }
2290 
2291  j++;
2292  }
2293 
2294  tol = resid * resid;
2295  for (i= 0; i<=m; i++)
2296  {
2297  delete v[i];
2298  }
2299  return 1;
2300 }
2301 
2303  const Operator *C_,
2304  const Operator *D_)
2305  : C(C_), D(D_), c_e(NULL), d_lo(NULL), d_hi(NULL), x_lo(NULL), x_hi(NULL),
2306  input_size(insize)
2307 {
2308  if (C) { MFEM_ASSERT(C->Width() == input_size, "Wrong width of C."); }
2309  if (D) { MFEM_ASSERT(D->Width() == input_size, "Wrong width of D."); }
2310 }
2311 
2313 {
2314  MFEM_ASSERT(C, "The C operator is unspecified -- can't set constraints.");
2315  MFEM_ASSERT(c.Size() == C->Height(), "Wrong size of the constraint.");
2316 
2317  c_e = &c;
2318 }
2319 
2321  const Vector &dh)
2322 {
2323  MFEM_ASSERT(D, "The D operator is unspecified -- can't set constraints.");
2324  MFEM_ASSERT(dl.Size() == D->Height() && dh.Size() == D->Height(),
2325  "Wrong size of the constraint.");
2326 
2327  d_lo = &dl; d_hi = &dh;
2328 }
2329 
2331 {
2332  MFEM_ASSERT(xl.Size() == input_size && xh.Size() == input_size,
2333  "Wrong size of the constraint.");
2334 
2335  x_lo = &xl; x_hi = &xh;
2336 }
2337 
2339 {
2340  int m = 0;
2341  if (C) { m += C->Height(); }
2342  if (D) { m += D->Height(); }
2343  return m;
2344 }
2345 
2347 {
2348  if (print_options.warnings)
2349  {
2350  MFEM_WARNING("Objective functional is ignored as SLBQP always minimizes"
2351  "the l2 norm of (x - x_target).");
2352  }
2353  MFEM_ASSERT(prob.GetC(), "Linear constraint is not set.");
2354  MFEM_ASSERT(prob.GetC()->Height() == 1, "Solver expects scalar constraint.");
2355 
2356  problem = &prob;
2357 }
2358 
2359 void SLBQPOptimizer::SetBounds(const Vector &lo_, const Vector &hi_)
2360 {
2361  lo.SetDataAndSize(lo_.GetData(), lo_.Size());
2362  hi.SetDataAndSize(hi_.GetData(), hi_.Size());
2363 }
2364 
2366 {
2367  w.SetDataAndSize(w_.GetData(), w_.Size());
2368  a = a_;
2369 }
2370 
2371 inline void SLBQPOptimizer::print_iteration(int it, double r, double l) const
2372 {
2374  {
2375  mfem::out << "SLBQP iteration " << it << ": residual = " << r
2376  << ", lambda = " << l << '\n';
2377  }
2378 }
2379 
2380 void SLBQPOptimizer::Mult(const Vector& xt, Vector& x) const
2381 {
2382  // Based on code provided by Denis Ridzal, dridzal@sandia.gov.
2383  // Algorithm adapted from Dai and Fletcher, "New Algorithms for
2384  // Singly Linearly Constrained Quadratic Programs Subject to Lower
2385  // and Upper Bounds", Numerical Analysis Report NA/216, 2003.
2386 
2387  // Set some algorithm-specific constants and temporaries.
2388  int nclip = 0;
2389  double l = 0;
2390  double llow = 0;
2391  double lupp = 0;
2392  double lnew = 0;
2393  double dl = 2;
2394  double r = 0;
2395  double rlow = 0;
2396  double rupp = 0;
2397  double s = 0;
2398 
2399  const double smin = 0.1;
2400 
2401  const double tol = max(abs_tol, rel_tol*a);
2402 
2403  // *** Start bracketing phase of SLBQP ***
2405  {
2406  mfem::out << "SLBQP bracketing phase" << '\n';
2407  }
2408 
2409  // Solve QP with fixed Lagrange multiplier
2410  r = initial_norm = solve(l,xt,x,nclip);
2411  print_iteration(nclip, r, l);
2412 
2413 
2414  // If x=xt was already within bounds and satisfies the linear
2415  // constraint, then we already have the solution.
2416  if (fabs(r) <= tol)
2417  {
2418  converged = true;
2419  goto slbqp_done;
2420  }
2421 
2422  if (r < 0)
2423  {
2424  llow = l; rlow = r; l = l + dl;
2425 
2426  // Solve QP with fixed Lagrange multiplier
2427  r = solve(l,xt,x,nclip);
2428  print_iteration(nclip, r, l);
2429 
2430  while ((r < 0) && (nclip < max_iter))
2431  {
2432  llow = l;
2433  s = rlow/r - 1.0;
2434  if (s < smin) { s = smin; }
2435  dl = dl + dl/s;
2436  l = l + dl;
2437 
2438  // Solve QP with fixed Lagrange multiplier
2439  r = solve(l,xt,x,nclip);
2440  print_iteration(nclip, r, l);
2441  }
2442 
2443  lupp = l; rupp = r;
2444  }
2445  else
2446  {
2447  lupp = l; rupp = r; l = l - dl;
2448 
2449  // Solve QP with fixed Lagrange multiplier
2450  r = solve(l,xt,x,nclip);
2451  print_iteration(nclip, r, l);
2452 
2453  while ((r > 0) && (nclip < max_iter))
2454  {
2455  lupp = l;
2456  s = rupp/r - 1.0;
2457  if (s < smin) { s = smin; }
2458  dl = dl + dl/s;
2459  l = l - dl;
2460 
2461  // Solve QP with fixed Lagrange multiplier
2462  r = solve(l,xt,x,nclip);
2463  print_iteration(nclip, r, l);
2464  }
2465 
2466  llow = l; rlow = r;
2467  }
2468 
2469  // *** Stop bracketing phase of SLBQP ***
2470 
2471 
2472  // *** Start secant phase of SLBQP ***
2474  {
2475  mfem::out << "SLBQP secant phase" << '\n';
2476  }
2477 
2478  s = 1.0 - rlow/rupp; dl = dl/s; l = lupp - dl;
2479 
2480  // Solve QP with fixed Lagrange multiplier
2481  r = solve(l,xt,x,nclip);
2482  print_iteration(nclip, r, l);
2483 
2484  while ( (fabs(r) > tol) && (nclip < max_iter) )
2485  {
2486  if (r > 0)
2487  {
2488  if (s <= 2.0)
2489  {
2490  lupp = l; rupp = r; s = 1.0 - rlow/rupp;
2491  dl = (lupp - llow)/s; l = lupp - dl;
2492  }
2493  else
2494  {
2495  s = rupp/r - 1.0;
2496  if (s < smin) { s = smin; }
2497  dl = (lupp - l)/s;
2498  lnew = 0.75*llow + 0.25*l;
2499  if (lnew < l-dl) { lnew = l-dl; }
2500  lupp = l; rupp = r; l = lnew;
2501  s = (lupp - llow)/(lupp - l);
2502  }
2503 
2504  }
2505  else
2506  {
2507  if (s >= 2.0)
2508  {
2509  llow = l; rlow = r; s = 1.0 - rlow/rupp;
2510  dl = (lupp - llow)/s; l = lupp - dl;
2511  }
2512  else
2513  {
2514  s = rlow/r - 1.0;
2515  if (s < smin) { s = smin; }
2516  dl = (l - llow)/s;
2517  lnew = 0.75*lupp + 0.25*l;
2518  if (lnew < l+dl) { lnew = l+dl; }
2519  llow = l; rlow = r; l = lnew;
2520  s = (lupp - llow)/(lupp - l);
2521  }
2522  }
2523 
2524  // Solve QP with fixed Lagrange multiplier
2525  r = solve(l,xt,x,nclip);
2526  print_iteration(nclip, r, l);
2527  }
2528 
2529  // *** Stop secant phase of SLBQP ***
2530  converged = (fabs(r) <= tol);
2531 
2532 slbqp_done:
2533 
2534  final_iter = nclip;
2535  final_norm = r;
2536 
2539  {
2540  mfem::out << "SLBQP: Number of iterations: " << final_iter << '\n'
2541  << " lambda = " << l << '\n'
2542  << " ||r|| = " << final_norm << '\n';
2543  }
2545  {
2546  mfem::out << "SLBQP: No convergence!" << '\n';
2547  }
2548 }
2549 
2550 struct WeightMinHeap
2551 {
2552  const std::vector<double> &w;
2553  std::vector<size_t> c;
2554  std::vector<int> loc;
2555 
2556  WeightMinHeap(const std::vector<double> &w_) : w(w_)
2557  {
2558  c.reserve(w.size());
2559  loc.resize(w.size());
2560  for (size_t i=0; i<w.size(); ++i) { push(i); }
2561  }
2562 
2563  size_t percolate_up(size_t pos, double val)
2564  {
2565  for (; pos > 0 && w[c[(pos-1)/2]] > val; pos = (pos-1)/2)
2566  {
2567  c[pos] = c[(pos-1)/2];
2568  loc[c[(pos-1)/2]] = pos;
2569  }
2570  return pos;
2571  }
2572 
2573  size_t percolate_down(size_t pos, double val)
2574  {
2575  while (2*pos+1 < c.size())
2576  {
2577  size_t left = 2*pos+1;
2578  size_t right = left+1;
2579  size_t tgt;
2580  if (right < c.size() && w[c[right]] < w[c[left]]) { tgt = right; }
2581  else { tgt = left; }
2582  if (w[c[tgt]] < val)
2583  {
2584  c[pos] = c[tgt];
2585  loc[c[tgt]] = pos;
2586  pos = tgt;
2587  }
2588  else
2589  {
2590  break;
2591  }
2592  }
2593  return pos;
2594  }
2595 
2596  void push(size_t i)
2597  {
2598  double val = w[i];
2599  c.push_back(0);
2600  size_t pos = c.size()-1;
2601  pos = percolate_up(pos, val);
2602  c[pos] = i;
2603  loc[i] = pos;
2604  }
2605 
2606  int pop()
2607  {
2608  size_t i = c[0];
2609  size_t j = c.back();
2610  c.pop_back();
2611  // Mark as removed
2612  loc[i] = -1;
2613  if (c.empty()) { return i; }
2614  double val = w[j];
2615  size_t pos = 0;
2616  pos = percolate_down(pos, val);
2617  c[pos] = j;
2618  loc[j] = pos;
2619  return i;
2620  }
2621 
2622  void update(size_t i)
2623  {
2624  size_t pos = loc[i];
2625  double val = w[i];
2626  pos = percolate_up(pos, val);
2627  pos = percolate_down(pos, val);
2628  c[pos] = i;
2629  loc[i] = pos;
2630  }
2631 
2632  bool picked(size_t i)
2633  {
2634  return loc[i] < 0;
2635  }
2636 };
2637 
2639 {
2640  int n = C.Width();
2641  // Scale rows by reciprocal of diagonal and take absolute value
2642  Vector D;
2643  C.GetDiag(D);
2644  int *I = C.GetI();
2645  int *J = C.GetJ();
2646  double *V = C.GetData();
2647  for (int i=0; i<n; ++i)
2648  {
2649  for (int j=I[i]; j<I[i+1]; ++j)
2650  {
2651  V[j] = abs(V[j]/D[i]);
2652  }
2653  }
2654 
2655  std::vector<double> w(n, 0.0);
2656  for (int k=0; k<n; ++k)
2657  {
2658  // Find all neighbors i of k
2659  for (int ii=I[k]; ii<I[k+1]; ++ii)
2660  {
2661  int i = J[ii];
2662  // Find value of (i,k)
2663  double C_ik = 0.0;
2664  for (int kk=I[i]; kk<I[i+1]; ++kk)
2665  {
2666  if (J[kk] == k)
2667  {
2668  C_ik = V[kk];
2669  break;
2670  }
2671  }
2672  for (int jj=I[k]; jj<I[k+1]; ++jj)
2673  {
2674  int j = J[jj];
2675  if (j == k) { continue; }
2676  double C_kj = V[jj];
2677  bool ij_exists = false;
2678  for (int jj2=I[i]; jj2<I[i+1]; ++jj2)
2679  {
2680  if (J[jj2] == j)
2681  {
2682  ij_exists = true;
2683  break;
2684  }
2685  }
2686  if (!ij_exists) { w[k] += pow(C_ik*C_kj,2); }
2687  }
2688  }
2689  w[k] = sqrt(w[k]);
2690  }
2691 
2692  WeightMinHeap w_heap(w);
2693 
2694  // Compute ordering
2695  p.SetSize(n);
2696  for (int ii=0; ii<n; ++ii)
2697  {
2698  int pi = w_heap.pop();
2699  p[ii] = pi;
2700  w[pi] = -1;
2701  for (int kk=I[pi]; kk<I[pi+1]; ++kk)
2702  {
2703  int k = J[kk];
2704  if (w_heap.picked(k)) { continue; }
2705  // Recompute weight
2706  w[k] = 0.0;
2707  // Find all neighbors i of k
2708  for (int ii2=I[k]; ii2<I[k+1]; ++ii2)
2709  {
2710  int i = J[ii2];
2711  if (w_heap.picked(i)) { continue; }
2712  // Find value of (i,k)
2713  double C_ik = 0.0;
2714  for (int kk2=I[i]; kk2<I[i+1]; ++kk2)
2715  {
2716  if (J[kk2] == k)
2717  {
2718  C_ik = V[kk2];
2719  break;
2720  }
2721  }
2722  for (int jj=I[k]; jj<I[k+1]; ++jj)
2723  {
2724  int j = J[jj];
2725  if (j == k || w_heap.picked(j)) { continue; }
2726  double C_kj = V[jj];
2727  bool ij_exists = false;
2728  for (int jj2=I[i]; jj2<I[i+1]; ++jj2)
2729  {
2730  if (J[jj2] == j)
2731  {
2732  ij_exists = true;
2733  break;
2734  }
2735  }
2736  if (!ij_exists) { w[k] += pow(C_ik*C_kj,2); }
2737  }
2738  }
2739  w[k] = sqrt(w[k]);
2740  w_heap.update(k);
2741  }
2742  }
2743 }
2744 
2745 BlockILU::BlockILU(int block_size_,
2746  Reordering reordering_,
2747  int k_fill_)
2748  : Solver(0),
2749  block_size(block_size_),
2750  k_fill(k_fill_),
2751  reordering(reordering_)
2752 { }
2753 
2755  int block_size_,
2756  Reordering reordering_,
2757  int k_fill_)
2758  : BlockILU(block_size_, reordering_, k_fill_)
2759 {
2760  SetOperator(op);
2761 }
2762 
2764 {
2765  const SparseMatrix *A = NULL;
2766 #ifdef MFEM_USE_MPI
2767  const HypreParMatrix *A_par = dynamic_cast<const HypreParMatrix *>(&op);
2768  SparseMatrix A_par_diag;
2769  if (A_par != NULL)
2770  {
2771  A_par->GetDiag(A_par_diag);
2772  A = &A_par_diag;
2773  }
2774 #endif
2775  if (A == NULL)
2776  {
2777  A = dynamic_cast<const SparseMatrix *>(&op);
2778  if (A == NULL)
2779  {
2780  MFEM_ABORT("BlockILU must be created with a SparseMatrix or HypreParMatrix");
2781  }
2782  }
2783  height = op.Height();
2784  width = op.Width();
2785  MFEM_ASSERT(A->Finalized(), "Matrix must be finalized.");
2786  CreateBlockPattern(*A);
2787  Factorize();
2788 }
2789 
2790 void BlockILU::CreateBlockPattern(const SparseMatrix &A)
2791 {
2792  MFEM_VERIFY(k_fill == 0, "Only block ILU(0) is currently supported.");
2793  if (A.Height() % block_size != 0)
2794  {
2795  MFEM_ABORT("BlockILU: block size must evenly divide the matrix size");
2796  }
2797 
2798  int nrows = A.Height();
2799  const int *I = A.GetI();
2800  const int *J = A.GetJ();
2801  const double *V = A.GetData();
2802  int nnz = 0;
2803  int nblockrows = nrows / block_size;
2804 
2805  std::vector<std::set<int>> unique_block_cols(nblockrows);
2806 
2807  for (int iblock = 0; iblock < nblockrows; ++iblock)
2808  {
2809  for (int bi = 0; bi < block_size; ++bi)
2810  {
2811  int i = iblock * block_size + bi;
2812  for (int k = I[i]; k < I[i + 1]; ++k)
2813  {
2814  unique_block_cols[iblock].insert(J[k] / block_size);
2815  }
2816  }
2817  nnz += unique_block_cols[iblock].size();
2818  }
2819 
2820  if (reordering != Reordering::NONE)
2821  {
2822  SparseMatrix C(nblockrows, nblockrows);
2823  for (int iblock = 0; iblock < nblockrows; ++iblock)
2824  {
2825  for (int jblock : unique_block_cols[iblock])
2826  {
2827  for (int bi = 0; bi < block_size; ++bi)
2828  {
2829  int i = iblock * block_size + bi;
2830  for (int k = I[i]; k < I[i + 1]; ++k)
2831  {
2832  int j = J[k];
2833  if (j >= jblock * block_size && j < (jblock + 1) * block_size)
2834  {
2835  C.Add(iblock, jblock, V[k]*V[k]);
2836  }
2837  }
2838  }
2839  }
2840  }
2841  C.Finalize(false);
2842  double *CV = C.GetData();
2843  for (int i=0; i<C.NumNonZeroElems(); ++i)
2844  {
2845  CV[i] = sqrt(CV[i]);
2846  }
2847 
2848  switch (reordering)
2849  {
2852  break;
2853  default:
2854  MFEM_ABORT("BlockILU: unknown reordering")
2855  }
2856  }
2857  else
2858  {
2859  // No reordering: permutation is identity
2860  P.SetSize(nblockrows);
2861  for (int i=0; i<nblockrows; ++i)
2862  {
2863  P[i] = i;
2864  }
2865  }
2866 
2867  // Compute inverse permutation
2868  Pinv.SetSize(nblockrows);
2869  for (int i=0; i<nblockrows; ++i)
2870  {
2871  Pinv[P[i]] = i;
2872  }
2873 
2874  // Permute columns
2875  std::vector<std::vector<int>> unique_block_cols_perminv(nblockrows);
2876  for (int i=0; i<nblockrows; ++i)
2877  {
2878  std::vector<int> &cols = unique_block_cols_perminv[i];
2879  for (int j : unique_block_cols[P[i]])
2880  {
2881  cols.push_back(Pinv[j]);
2882  }
2883  std::sort(cols.begin(), cols.end());
2884  }
2885 
2886  ID.SetSize(nblockrows);
2887  IB.SetSize(nblockrows + 1);
2888  IB[0] = 0;
2889  JB.SetSize(nnz);
2890  AB.SetSize(block_size, block_size, nnz);
2891  DB.SetSize(block_size, block_size, nblockrows);
2892  AB = 0.0;
2893  DB = 0.0;
2894  ipiv.SetSize(block_size*nblockrows);
2895  int counter = 0;
2896 
2897  for (int iblock = 0; iblock < nblockrows; ++iblock)
2898  {
2899  int iblock_perm = P[iblock];
2900  for (int jblock : unique_block_cols_perminv[iblock])
2901  {
2902  int jblock_perm = P[jblock];
2903  if (iblock == jblock)
2904  {
2905  ID[iblock] = counter;
2906  }
2907  JB[counter] = jblock;
2908  for (int bi = 0; bi < block_size; ++bi)
2909  {
2910  int i = iblock_perm*block_size + bi;
2911  for (int k = I[i]; k < I[i + 1]; ++k)
2912  {
2913  int j = J[k];
2914  if (j >= jblock_perm*block_size && j < (jblock_perm + 1)*block_size)
2915  {
2916  int bj = j - jblock_perm*block_size;
2917  double val = V[k];
2918  AB(bi, bj, counter) = val;
2919  // Extract the diagonal
2920  if (iblock == jblock)
2921  {
2922  DB(bi, bj, iblock) = val;
2923  }
2924  }
2925  }
2926  }
2927  ++counter;
2928  }
2929  IB[iblock + 1] = counter;
2930  }
2931 }
2932 
2933 void BlockILU::Factorize()
2934 {
2935  int nblockrows = Height()/block_size;
2936 
2937  // Precompute LU factorization of diagonal blocks
2938  for (int i=0; i<nblockrows; ++i)
2939  {
2940  LUFactors factorization(DB.GetData(i), &ipiv[i*block_size]);
2941  factorization.Factor(block_size);
2942  }
2943 
2944  // Note: we use UseExternalData to extract submatrices from the tensor AB
2945  // instead of the DenseTensor call operator, because the call operator does
2946  // not allow for two simultaneous submatrix views into the same tensor
2947  DenseMatrix A_ik, A_ij, A_kj;
2948  // Loop over block rows (starting with second block row)
2949  for (int i=1; i<nblockrows; ++i)
2950  {
2951  // Find all nonzeros to the left of the diagonal in row i
2952  for (int kk=IB[i]; kk<IB[i+1]; ++kk)
2953  {
2954  int k = JB[kk];
2955  // Make sure we're still to the left of the diagonal
2956  if (k == i) { break; }
2957  if (k > i)
2958  {
2959  MFEM_ABORT("Matrix must be sorted with nonzero diagonal");
2960  }
2961  LUFactors A_kk_inv(DB.GetData(k), &ipiv[k*block_size]);
2962  A_ik.UseExternalData(&AB(0,0,kk), block_size, block_size);
2963  // A_ik = A_ik * A_kk^{-1}
2964  A_kk_inv.RightSolve(block_size, block_size, A_ik.GetData());
2965  // Modify everything to the right of k in row i
2966  for (int jj=kk+1; jj<IB[i+1]; ++jj)
2967  {
2968  int j = JB[jj];
2969  if (j <= k) { continue; } // Superfluous because JB is sorted?
2970  A_ij.UseExternalData(&AB(0,0,jj), block_size, block_size);
2971  for (int ll=IB[k]; ll<IB[k+1]; ++ll)
2972  {
2973  int l = JB[ll];
2974  if (l == j)
2975  {
2976  A_kj.UseExternalData(&AB(0,0,ll), block_size, block_size);
2977  // A_ij = A_ij - A_ik*A_kj;
2978  AddMult_a(-1.0, A_ik, A_kj, A_ij);
2979  // If we need to, update diagonal factorization
2980  if (j == i)
2981  {
2982  DB(i) = A_ij;
2983  LUFactors factorization(DB.GetData(i), &ipiv[i*block_size]);
2984  factorization.Factor(block_size);
2985  }
2986  break;
2987  }
2988  }
2989  }
2990  }
2991  }
2992 }
2993 
2994 void BlockILU::Mult(const Vector &b, Vector &x) const
2995 {
2996  MFEM_ASSERT(height > 0, "BlockILU(0) preconditioner is not constructed");
2997  int nblockrows = Height()/block_size;
2998  y.SetSize(Height());
2999 
3000  DenseMatrix B;
3001  Vector yi, yj, xi, xj;
3002  Vector tmp(block_size);
3003  // Forward substitute to solve Ly = b
3004  // Implicitly, L has identity on the diagonal
3005  y = 0.0;
3006  for (int i=0; i<nblockrows; ++i)
3007  {
3008  yi.SetDataAndSize(&y[i*block_size], block_size);
3009  for (int ib=0; ib<block_size; ++ib)
3010  {
3011  yi[ib] = b[ib + P[i]*block_size];
3012  }
3013  for (int k=IB[i]; k<ID[i]; ++k)
3014  {
3015  int j = JB[k];
3016  const DenseMatrix &L_ij = AB(k);
3017  yj.SetDataAndSize(&y[j*block_size], block_size);
3018  // y_i = y_i - L_ij*y_j
3019  L_ij.AddMult_a(-1.0, yj, yi);
3020  }
3021  }
3022  // Backward substitution to solve Ux = y
3023  for (int i=nblockrows-1; i >= 0; --i)
3024  {
3025  xi.SetDataAndSize(&x[P[i]*block_size], block_size);
3026  for (int ib=0; ib<block_size; ++ib)
3027  {
3028  xi[ib] = y[ib + i*block_size];
3029  }
3030  for (int k=ID[i]+1; k<IB[i+1]; ++k)
3031  {
3032  int j = JB[k];
3033  const DenseMatrix &U_ij = AB(k);
3034  xj.SetDataAndSize(&x[P[j]*block_size], block_size);
3035  // x_i = x_i - U_ij*x_j
3036  U_ij.AddMult_a(-1.0, xj, xi);
3037  }
3038  LUFactors A_ii_inv(&DB(0,0,i), &ipiv[i*block_size]);
3039  // x_i = D_ii^{-1} x_i
3040  A_ii_inv.Solve(block_size, 1, xi.GetData());
3041  }
3042 }
3043 
3044 
3046  int it, double norm, const Vector &r, bool final)
3047 {
3048  if (!ess_dofs_list) { return; }
3049 
3050  double bc_norm_squared = 0.0;
3051  r.HostRead();
3053  for (int i = 0; i < ess_dofs_list->Size(); i++)
3054  {
3055  const double r_entry = r((*ess_dofs_list)[i]);
3056  bc_norm_squared += r_entry*r_entry;
3057  }
3058  bool print = true;
3059 #ifdef MFEM_USE_MPI
3060  MPI_Comm comm = iter_solver->GetComm();
3061  if (comm != MPI_COMM_NULL)
3062  {
3063  double glob_bc_norm_squared = 0.0;
3064  MPI_Reduce(&bc_norm_squared, &glob_bc_norm_squared, 1, MPI_DOUBLE,
3065  MPI_SUM, 0, comm);
3066  bc_norm_squared = glob_bc_norm_squared;
3067  int rank;
3068  MPI_Comm_rank(comm, &rank);
3069  print = (rank == 0);
3070  }
3071 #endif
3072  if ((it == 0 || final || bc_norm_squared > 0.0) && print)
3073  {
3074  mfem::out << " ResidualBCMonitor : b.c. residual norm = "
3075  << sqrt(bc_norm_squared) << endl;
3076  }
3077 }
3078 
3079 
3080 #ifdef MFEM_USE_SUITESPARSE
3081 
3083 {
3084  mat = NULL;
3085  Numeric = NULL;
3086  AI = AJ = NULL;
3087  if (!use_long_ints)
3088  {
3089  umfpack_di_defaults(Control);
3090  }
3091  else
3092  {
3093  umfpack_dl_defaults(Control);
3094  }
3095 }
3096 
3098 {
3099  void *Symbolic;
3100 
3101  if (Numeric)
3102  {
3103  if (!use_long_ints)
3104  {
3105  umfpack_di_free_numeric(&Numeric);
3106  }
3107  else
3108  {
3109  umfpack_dl_free_numeric(&Numeric);
3110  }
3111  }
3112 
3113  mat = const_cast<SparseMatrix *>(dynamic_cast<const SparseMatrix *>(&op));
3114  MFEM_VERIFY(mat, "not a SparseMatrix");
3115 
3116  // UMFPack requires that the column-indices in mat corresponding to each
3117  // row be sorted.
3118  // Generally, this will modify the ordering of the entries of mat.
3120 
3121  height = mat->Height();
3122  width = mat->Width();
3123  MFEM_VERIFY(width == height, "not a square matrix");
3124 
3125  const int * Ap = mat->HostReadI();
3126  const int * Ai = mat->HostReadJ();
3127  const double * Ax = mat->HostReadData();
3128 
3129  if (!use_long_ints)
3130  {
3131  int status = umfpack_di_symbolic(width, width, Ap, Ai, Ax, &Symbolic,
3132  Control, Info);
3133  if (status < 0)
3134  {
3135  umfpack_di_report_info(Control, Info);
3136  umfpack_di_report_status(Control, status);
3137  mfem_error("UMFPackSolver::SetOperator :"
3138  " umfpack_di_symbolic() failed!");
3139  }
3140 
3141  status = umfpack_di_numeric(Ap, Ai, Ax, Symbolic, &Numeric,
3142  Control, Info);
3143  if (status < 0)
3144  {
3145  umfpack_di_report_info(Control, Info);
3146  umfpack_di_report_status(Control, status);
3147  mfem_error("UMFPackSolver::SetOperator :"
3148  " umfpack_di_numeric() failed!");
3149  }
3150  umfpack_di_free_symbolic(&Symbolic);
3151  }
3152  else
3153  {
3154  SuiteSparse_long status;
3155 
3156  delete [] AJ;
3157  delete [] AI;
3158  AI = new SuiteSparse_long[width + 1];
3159  AJ = new SuiteSparse_long[Ap[width]];
3160  for (int i = 0; i <= width; i++)
3161  {
3162  AI[i] = (SuiteSparse_long)(Ap[i]);
3163  }
3164  for (int i = 0; i < Ap[width]; i++)
3165  {
3166  AJ[i] = (SuiteSparse_long)(Ai[i]);
3167  }
3168 
3169  status = umfpack_dl_symbolic(width, width, AI, AJ, Ax, &Symbolic,
3170  Control, Info);
3171  if (status < 0)
3172  {
3173  umfpack_dl_report_info(Control, Info);
3174  umfpack_dl_report_status(Control, status);
3175  mfem_error("UMFPackSolver::SetOperator :"
3176  " umfpack_dl_symbolic() failed!");
3177  }
3178 
3179  status = umfpack_dl_numeric(AI, AJ, Ax, Symbolic, &Numeric,
3180  Control, Info);
3181  if (status < 0)
3182  {
3183  umfpack_dl_report_info(Control, Info);
3184  umfpack_dl_report_status(Control, status);
3185  mfem_error("UMFPackSolver::SetOperator :"
3186  " umfpack_dl_numeric() failed!");
3187  }
3188  umfpack_dl_free_symbolic(&Symbolic);
3189  }
3190 }
3191 
3192 void UMFPackSolver::Mult(const Vector &b, Vector &x) const
3193 {
3194  if (mat == NULL)
3195  mfem_error("UMFPackSolver::Mult : matrix is not set!"
3196  " Call SetOperator first!");
3197  b.HostRead();
3198  x.HostReadWrite();
3199  if (!use_long_ints)
3200  {
3201  int status =
3202  umfpack_di_solve(UMFPACK_At, mat->HostReadI(), mat->HostReadJ(),
3203  mat->HostReadData(), x.HostWrite(), b.HostRead(),
3204  Numeric, Control, Info);
3205  umfpack_di_report_info(Control, Info);
3206  if (status < 0)
3207  {
3208  umfpack_di_report_status(Control, status);
3209  mfem_error("UMFPackSolver::Mult : umfpack_di_solve() failed!");
3210  }
3211  }
3212  else
3213  {
3214  SuiteSparse_long status =
3215  umfpack_dl_solve(UMFPACK_At, AI, AJ, mat->HostReadData(),
3216  x.HostWrite(), b.HostRead(), Numeric, Control,
3217  Info);
3218  umfpack_dl_report_info(Control, Info);
3219  if (status < 0)
3220  {
3221  umfpack_dl_report_status(Control, status);
3222  mfem_error("UMFPackSolver::Mult : umfpack_dl_solve() failed!");
3223  }
3224  }
3225 }
3226 
3228 {
3229  if (mat == NULL)
3230  mfem_error("UMFPackSolver::MultTranspose : matrix is not set!"
3231  " Call SetOperator first!");
3232  b.HostRead();
3233  x.HostReadWrite();
3234  if (!use_long_ints)
3235  {
3236  int status =
3237  umfpack_di_solve(UMFPACK_A, mat->HostReadI(), mat->HostReadJ(),
3238  mat->HostReadData(), x.HostWrite(), b.HostRead(),
3239  Numeric, Control, Info);
3240  umfpack_di_report_info(Control, Info);
3241  if (status < 0)
3242  {
3243  umfpack_di_report_status(Control, status);
3244  mfem_error("UMFPackSolver::MultTranspose :"
3245  " umfpack_di_solve() failed!");
3246  }
3247  }
3248  else
3249  {
3250  SuiteSparse_long status =
3251  umfpack_dl_solve(UMFPACK_A, AI, AJ, mat->HostReadData(),
3252  x.HostWrite(), b.HostRead(), Numeric, Control,
3253  Info);
3254  umfpack_dl_report_info(Control, Info);
3255  if (status < 0)
3256  {
3257  umfpack_dl_report_status(Control, status);
3258  mfem_error("UMFPackSolver::MultTranspose :"
3259  " umfpack_dl_solve() failed!");
3260  }
3261  }
3262 }
3263 
3265 {
3266  delete [] AJ;
3267  delete [] AI;
3268  if (Numeric)
3269  {
3270  if (!use_long_ints)
3271  {
3272  umfpack_di_free_numeric(&Numeric);
3273  }
3274  else
3275  {
3276  umfpack_dl_free_numeric(&Numeric);
3277  }
3278  }
3279 }
3280 
3282 {
3283  klu_defaults(&Common);
3284 }
3285 
3287 {
3288  if (Numeric)
3289  {
3290  MFEM_ASSERT(Symbolic != 0,
3291  "Had Numeric pointer in KLU, but not Symbolic");
3292  klu_free_symbolic(&Symbolic, &Common);
3293  Symbolic = 0;
3294  klu_free_numeric(&Numeric, &Common);
3295  Numeric = 0;
3296  }
3297 
3298  mat = const_cast<SparseMatrix *>(dynamic_cast<const SparseMatrix *>(&op));
3299  MFEM_VERIFY(mat != NULL, "not a SparseMatrix");
3300 
3301  // KLU requires that the column-indices in mat corresponding to each row be
3302  // sorted. Generally, this will modify the ordering of the entries of mat.
3304 
3305  height = mat->Height();
3306  width = mat->Width();
3307  MFEM_VERIFY(width == height, "not a square matrix");
3308 
3309  int * Ap = mat->GetI();
3310  int * Ai = mat->GetJ();
3311  double * Ax = mat->GetData();
3312 
3313  Symbolic = klu_analyze( height, Ap, Ai, &Common);
3314  Numeric = klu_factor(Ap, Ai, Ax, Symbolic, &Common);
3315 }
3316 
3317 void KLUSolver::Mult(const Vector &b, Vector &x) const
3318 {
3319  MFEM_VERIFY(mat != NULL,
3320  "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
3321 
3322  int n = mat->Height();
3323  int numRhs = 1;
3324  // Copy B into X, so we can pass it in and overwrite it.
3325  x = b;
3326  // Solve the transpose, since KLU thinks the matrix is compressed column
3327  // format.
3328  klu_tsolve( Symbolic, Numeric, n, numRhs, x.GetData(), &Common);
3329 }
3330 
3331 void KLUSolver::MultTranspose(const Vector &b, Vector &x) const
3332 {
3333  MFEM_VERIFY(mat != NULL,
3334  "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
3335 
3336  int n = mat->Height();
3337  int numRhs = 1;
3338  // Copy B into X, so we can pass it in and overwrite it.
3339  x = b;
3340  // Solve the regular matrix, not the transpose, since KLU thinks the matrix
3341  // is compressed column format.
3342  klu_solve( Symbolic, Numeric, n, numRhs, x.GetData(), &Common);
3343 }
3344 
3346 {
3347  klu_free_symbolic (&Symbolic, &Common) ;
3348  klu_free_numeric (&Numeric, &Common) ;
3349  Symbolic = 0;
3350  Numeric = 0;
3351 }
3352 
3353 #endif // MFEM_USE_SUITESPARSE
3354 
3356  const SparseMatrix &block_dof_)
3357  : Solver(A.NumRows()), block_dof(const_cast<SparseMatrix&>(block_dof_)),
3358  block_solvers(new DenseMatrixInverse[block_dof.NumRows()])
3359 {
3360  DenseMatrix sub_A;
3361  for (int i = 0; i < block_dof.NumRows(); ++i)
3362  {
3363  local_dofs.MakeRef(block_dof.GetRowColumns(i), block_dof.RowSize(i));
3364  sub_A.SetSize(local_dofs.Size());
3365  A.GetSubMatrix(local_dofs, local_dofs, sub_A);
3366  block_solvers[i].SetOperator(sub_A);
3367  }
3368 }
3369 
3370 void DirectSubBlockSolver::Mult(const Vector &x, Vector &y) const
3371 {
3372  y.SetSize(x.Size());
3373  y = 0.0;
3374 
3375  for (int i = 0; i < block_dof.NumRows(); ++i)
3376  {
3377  local_dofs.MakeRef(block_dof.GetRowColumns(i), block_dof.RowSize(i));
3378  x.GetSubVector(local_dofs, sub_rhs);
3379  sub_sol.SetSize(local_dofs.Size());
3380  block_solvers[i].Mult(sub_rhs, sub_sol);
3381  y.AddElementVector(local_dofs, sub_sol);
3382  }
3383 }
3384 
3385 void ProductSolver::Mult(const Vector & x, Vector & y) const
3386 {
3387  y.SetSize(x.Size());
3388  y = 0.0;
3389  S0->Mult(x, y);
3390 
3391  Vector z(x.Size());
3392  z = 0.0;
3393  A->Mult(y, z);
3394  add(-1.0, z, 1.0, x, z); // z = (I - A * S0) x
3395 
3396  Vector S1z(x.Size());
3397  S1z = 0.0;
3398  S1->Mult(z, S1z);
3399  y += S1z;
3400 }
3401 
3402 void ProductSolver::MultTranspose(const Vector & x, Vector & y) const
3403 {
3404  y.SetSize(x.Size());
3405  y = 0.0;
3406  S1->MultTranspose(x, y);
3407 
3408  Vector z(x.Size());
3409  z = 0.0;
3410  A->MultTranspose(y, z);
3411  add(-1.0, z, 1.0, x, z); // z = (I - A^T * S1^T) x
3412 
3413  Vector S0Tz(x.Size());
3414  S0Tz = 0.0;
3415  S0->MultTranspose(z, S0Tz);
3416  y += S0Tz;
3417 }
3418 
3420  : Solver(0, false), global_size(-1)
3421 #ifdef MFEM_USE_MPI
3422  , parallel(false)
3423 #endif
3424 { }
3425 
3426 #ifdef MFEM_USE_MPI
3427 OrthoSolver::OrthoSolver(MPI_Comm mycomm_)
3428  : Solver(0, false), mycomm(mycomm_), global_size(-1), parallel(true) { }
3429 #endif
3430 
3432 {
3433  solver = &s;
3434  height = s.Height();
3435  width = s.Width();
3436  MFEM_VERIFY(height == width, "Solver must be a square Operator!");
3437  global_size = -1; // lazy evaluated
3438 }
3439 
3441 {
3442  MFEM_VERIFY(solver, "Solver hasn't been set, call SetSolver() first.");
3443  solver->SetOperator(op);
3444  height = solver->Height();
3445  width = solver->Width();
3446  MFEM_VERIFY(height == width, "Solver must be a square Operator!");
3447  global_size = -1; // lazy evaluated
3448 }
3449 
3450 void OrthoSolver::Mult(const Vector &b, Vector &x) const
3451 {
3452  MFEM_VERIFY(solver, "Solver hasn't been set, call SetSolver() first.");
3453  MFEM_VERIFY(height == solver->Height(),
3454  "solver was modified externally! call SetSolver() again!");
3455  MFEM_VERIFY(height == b.Size(), "incompatible input Vector size!");
3456  MFEM_VERIFY(height == x.Size(), "incompatible output Vector size!");
3457 
3458  // Orthogonalize input
3459  Orthogonalize(b, b_ortho);
3460 
3461  // Propagate iterative_mode to the solver:
3462  solver->iterative_mode = iterative_mode;
3463 
3464  // Apply the Solver
3465  solver->Mult(b_ortho, x);
3466 
3467  // Orthogonalize output
3468  Orthogonalize(x, x);
3469 }
3470 
3471 void OrthoSolver::Orthogonalize(const Vector &v, Vector &v_ortho) const
3472 {
3473  if (global_size == -1)
3474  {
3475  global_size = height;
3476 #ifdef MFEM_USE_MPI
3477  if (parallel)
3478  {
3479  MPI_Allreduce(MPI_IN_PLACE, &global_size, 1, HYPRE_MPI_BIG_INT,
3480  MPI_SUM, mycomm);
3481  }
3482 #endif
3483  }
3484 
3485  // TODO: GPU/device implementation
3486 
3487  double global_sum = v.Sum();
3488 
3489 #ifdef MFEM_USE_MPI
3490  if (parallel)
3491  {
3492  MPI_Allreduce(MPI_IN_PLACE, &global_sum, 1, MPI_DOUBLE, MPI_SUM, mycomm);
3493  }
3494 #endif
3495 
3496  double ratio = global_sum / static_cast<double>(global_size);
3497  v_ortho.SetSize(v.Size());
3498  v.HostRead();
3499  v_ortho.HostWrite();
3500  for (int i = 0; i < v_ortho.Size(); ++i)
3501  {
3502  v_ortho(i) = v(i) - ratio;
3503  }
3504 }
3505 
3506 #ifdef MFEM_USE_MPI
3508  HypreParMatrix *aux_map,
3509  bool op_is_symmetric,
3510  bool own_aux_map)
3511  : Solver(op.NumRows()), aux_map_(aux_map, own_aux_map)
3512 {
3513  aux_system_.Reset(RAP(&op, aux_map));
3514  aux_system_.As<HypreParMatrix>()->EliminateZeroRows();
3515  aux_smoother_.Reset(new HypreSmoother(*aux_system_.As<HypreParMatrix>()));
3516  aux_smoother_.As<HypreSmoother>()->SetOperatorSymmetry(op_is_symmetric);
3517 }
3518 
3519 void AuxSpaceSmoother::Mult(const Vector &x, Vector &y, bool transpose) const
3520 {
3521  Vector aux_rhs(aux_map_->NumCols());
3522  aux_map_->MultTranspose(x, aux_rhs);
3523 
3524  Vector aux_sol(aux_rhs.Size());
3525  if (transpose)
3526  {
3527  aux_smoother_->MultTranspose(aux_rhs, aux_sol);
3528  }
3529  else
3530  {
3531  aux_smoother_->Mult(aux_rhs, aux_sol);
3532  }
3533 
3534  y.SetSize(aux_map_->NumRows());
3535  aux_map_->Mult(aux_sol, y);
3536 }
3537 #endif // MFEM_USE_MPI
3538 
3539 }
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:307
const Array< int > * ess_dofs_list
Not owned.
Definition: solvers.hpp:1039
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: solvers.cpp:946
PowerMethod helper class to estimate the largest eigenvalue of an operator using the iterative power ...
Definition: operator.hpp:986
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:93
const double * HostReadData() const
Definition: sparsemat.hpp:277
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:733
Conjugate gradient method.
Definition: solvers.hpp:493
double Info[UMFPACK_INFO]
Definition: solvers.hpp:1065
Chebyshev accelerated smoothing with given vector, no matrix necessary.
Definition: solvers.hpp:379
SparseMatrix * mat
Definition: solvers.hpp:1057
const Operator * oper
Definition: solvers.hpp:120
Array< Vector * > ykArray
Definition: solvers.hpp:730
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
Definition: operator.hpp:122
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:669
SuiteSparse_long * AI
Definition: solvers.hpp:1059
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:512
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:452
void MinimumDiscardedFillOrdering(SparseMatrix &C, Array< int > &p)
Definition: solvers.cpp:2638
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:391
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:344
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:117
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
Definition: solvers.hpp:91
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:300
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
Definition: solvers.hpp:704
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:2005
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
Definition: solvers.hpp:94
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1162
OperatorJacobiSmoother(const double damping=1.0)
Default constructor: the diagonal will be computed by subsequent calls to SetOperator() using the Ope...
Definition: solvers.cpp:200
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:227
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1612
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:712
virtual void Solve(int m, int n, double *X) const
Definition: densemat.cpp:3420
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:974
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:460
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: operator.hpp:86
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:525
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:448
int * GetI()
Return the array I.
Definition: sparsemat.hpp:222
bool errors
If a fatal problem has been detected the failure will be reported to mfem::err.
Definition: solvers.hpp:82
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
Definition: operator.hpp:69
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3317
void SetOperator(const Operator &op)
Definition: solvers.cpp:2763
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:686
STL namespace.
MINRES method.
Definition: solvers.hpp:603
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:548
Reordering
The reordering method used by the BlockILU factorization.
Definition: solvers.hpp:951
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2886
void Mult(const Vector &b, Vector &x) const
Solve the system LUx = b, where L and U are the block ILU factors.
Definition: solvers.cpp:2994
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
Definition: solvers.hpp:130
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:457
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:314
void Monitor(int it, double norm, const Vector &r, const Vector &x, bool final=false) const
Definition: solvers.cpp:190
const Vector * d_hi
Definition: solvers.hpp:818
double GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
Definition: solvers.hpp:265
virtual void SetOperator(const Operator &op)
Set the Operator that is the OrthoSolver is to invert (approximately).
Definition: solvers.cpp:3440
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:232
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:3370
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:1568
bool IsFinite(const double &val)
Definition: vector.hpp:481
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:179
const Vector * x_lo
Definition: solvers.hpp:818
double Norm(const Vector &x) const
Definition: solvers.hpp:172
Operator * grad
Definition: solvers.hpp:645
const Vector * x_hi
Definition: solvers.hpp:818
const class IterativeSolver * iter_solver
The last IterativeSolver to which this monitor was attached.
Definition: solvers.hpp:40
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2321
General product operator: x -> (A*B)(x) = A(B(x)).
Definition: operator.hpp:774
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:616
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:3331
int GetNumConstraints() const
Definition: solvers.cpp:2338
Data type sparse matrix.
Definition: sparsemat.hpp:50
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Definition: solvers.hpp:887
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:1778
virtual ~UMFPackSolver()
Definition: solvers.cpp:3264
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:456
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:581
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition: solvers.hpp:250
double b
Definition: lissajous.cpp:42
void Update(Vector &x, int k, DenseMatrix &h, Vector &s, Array< Vector *> &v)
Definition: solvers.cpp:953
const int * HostReadJ() const
Definition: sparsemat.hpp:261
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double delta
Definition: lissajous.cpp:43
OperatorChebyshevSmoother(const Operator &oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, double max_eig_estimate)
Definition: solvers.cpp:328
void SetLinearConstraint(const Vector &w_, double a_)
Definition: solvers.cpp:2365
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:912
DirectSubBlockSolver(const SparseMatrix &A, const SparseMatrix &block_dof)
Definition: solvers.cpp:3355
SuiteSparse_long * AJ
Definition: solvers.hpp:1059
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition: solvers.hpp:538
void print_iteration(int it, double r, double l) const
Definition: solvers.cpp:2371
void AddMult_a(double a, const Vector &x, Vector &y) const
y += a * A.x
Definition: densemat.cpp:299
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:892
double * GetData(int k)
Definition: densemat.hpp:1083
void Setup(const Vector &diag)
Definition: solvers.cpp:277
Parallel smoothers in hypre.
Definition: hypre.hpp:971
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:907
void SetSize(int i, int j, int k, MemoryType mt_=MemoryType::PRESERVE)
Definition: densemat.hpp:1030
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:3286
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1591
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:1481
Array< Vector * > skArray
Definition: solvers.hpp:730
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Definition: densemat.cpp:3252
prob_type prob
Definition: ex25.cpp:153
klu_symbolic * Symbolic
Definition: solvers.hpp:1098
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:311
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:640
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:552
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
Definition: solvers.cpp:2745
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
Stationary linear iteration: x <- x + B (b - A x)
Definition: solvers.hpp:461
double omega
Definition: ex25.cpp:141
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
Definition: solvers.cpp:2320
void SetEqualityConstraint(const Vector &c)
Definition: solvers.cpp:2312
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:590
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:908
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1817
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:635
void SetAbsTol(double atol)
Definition: solvers.hpp:200
PrintLevel print_options
Output behavior for the iterative solver.
Definition: solvers.hpp:137
void SetRelTol(double rtol)
Definition: solvers.hpp:199
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:1064
Abstract base class for iterative solver.
Definition: solvers.hpp:66
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:3227
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:461
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
void SetAdaptiveLinRtol(const int type=2, const double rtol0=0.5, const double rtol_max=0.9, const double alpha=0.5 *(1.0+sqrt(5.0)), const double gamma=1.0)
Enable adaptive linear solver relative tolerance algorithm.
Definition: solvers.cpp:1925
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:684
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
Definition: solvers.hpp:252
int max_iter
Limit for the number of iterations the solver is allowed to do.
Definition: solvers.hpp:151
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3192
GMRES method.
Definition: solvers.hpp:525
void AdaptiveLinRtolPostSolve(const Vector &x, const Vector &b, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked after the linear solve.
Definition: solvers.cpp:1985
static int GuessLegacyPrintLevel(PrintLevel)
Use some heuristics to guess a legacy print level corresponding to the given PrintLevel.
Definition: solvers.cpp:149
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int NumCols() const
Get the number of columns (size of input) of the Operator. Synonym with Width().
Definition: operator.hpp:75
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:156
const int * HostReadI() const
Definition: sparsemat.hpp:245
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:268
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.cpp:2346
double a
Definition: lissajous.cpp:41
void UpdateVectors()
Definition: solvers.cpp:519
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:414
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:464
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:253
klu_common Common
Definition: solvers.hpp:1119
double rel_tol
Relative tolerance.
Definition: solvers.hpp:154
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
Definition: solvers.cpp:2302
void SetSolutionBounds(const Vector &xl, const Vector &xh)
Definition: solvers.cpp:2330
virtual void Mult(const Vector &xt, Vector &x) const
Definition: solvers.cpp:2380
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
AuxSpaceSmoother(const HypreParMatrix &op, HypreParMatrix *aux_map, bool op_is_symmetric=true, bool own_aux_map=false)
Definition: solvers.cpp:3507
void SetOperator(const Operator &op)
Recompute the diagonal using the method AssembleDiagonal of the given new Operator, op.
Definition: solvers.cpp:241
void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: solvers.cpp:924
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
Definition: mem_manager.cpp:55
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Definition: operator.hpp:131
PrintLevel FromLegacyPrintLevel(int)
Convert a legacy print level integer to a PrintLevel object.
Definition: solvers.cpp:114
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:586
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:804
double Dot(const Vector &x, const Vector &y) const
Definition: solvers.cpp:55
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
Definition: solvers.hpp:292
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:869
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Definition: solvers.cpp:3431
void Mult(const Vector &b, Vector &x) const
Perform the action of the OrthoSolver: P * solver * P where P is the projection to the subspace of ve...
Definition: solvers.cpp:3450
void MonitorResidual(int it, double norm, const Vector &r, bool final) override
Monitor the residual vector r.
Definition: solvers.cpp:3045
const Vector * d_lo
Definition: solvers.hpp:818
virtual void MonitorSolution(int it, double norm, const Vector &x, bool final)
Monitor the solution vector x.
Definition: solvers.hpp:54
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
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:2147
virtual ~KLUSolver()
Definition: solvers.cpp:3345
RefCoord s[3]
void UpdateVectors()
Definition: solvers.cpp:703
BiCGSTAB method.
Definition: solvers.hpp:572
bool warnings
If a non-fatal problem has been detected some context-specific information will be reported to mfem::...
Definition: solvers.hpp:85
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition: solvers.hpp:88
SparseMatrix * mat
Definition: solvers.hpp:1097
Base class for solvers.
Definition: operator.hpp:682
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
const OptimizationProblem * problem
Definition: solvers.hpp:851
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1376
double EstimateLargestEigenvalue(Operator &opr, Vector &v0, int numSteps=10, double tolerance=1e-8, int seed=12345)
Returns an estimate of the largest eigenvalue of the operator opr using the iterative power method...
Definition: operator.cpp:695
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
Definition: solvers.hpp:699
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
void SetBounds(const Vector &lo_, const Vector &hi_)
Definition: solvers.cpp:2359
void AddMult_a(double alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
Definition: densemat.cpp:2414
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:3385
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:3402
void AdaptiveLinRtolPreSolve(const Vector &x, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked before the linear solve.
Definition: solvers.cpp:1938
const Operator * D
Definition: solvers.hpp:817
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:479
const Vector * c_e
Definition: solvers.hpp:818
IterativeSolverMonitor * monitor
Definition: solvers.hpp:122
klu_numeric * Numeric
Definition: solvers.hpp:1099
virtual void MonitorResidual(int it, double norm, const Vector &r, bool final)
Monitor the residual vector r.
Definition: solvers.hpp:48
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:468
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:475
Settings for the output behavior of the IterativeSolver.
Definition: solvers.hpp:78
const Operator * C
Not owned, some can remain unused (NULL).
Definition: solvers.hpp:817
double abs_tol
Absolute tolerance.
Definition: solvers.hpp:157
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
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:671
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:1338
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1806
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:3097