MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
hypre.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "linalg.hpp"
17 #include "../fem/fem.hpp"
18 
19 #include <fstream>
20 #include <iomanip>
21 #include <cmath>
22 #include <cstdlib>
23 
24 using namespace std;
25 
26 namespace mfem
27 {
28 
29 template<typename TargetT, typename SourceT>
30 static TargetT *DuplicateAs(const SourceT *array, int size,
31  bool cplusplus = true)
32 {
33  TargetT *target_array = cplusplus ? new TargetT[size]
34  /* */ : hypre_TAlloc(TargetT, size);
35  for (int i = 0; i < size; i++)
36  {
37  target_array[i] = array[i];
38  }
39  return target_array;
40 }
41 
42 inline void HypreParVector::_SetDataAndSize_()
43 {
44  SetDataAndSize(hypre_VectorData(hypre_ParVectorLocalVector(x)),
46  hypre_VectorSize(hypre_ParVectorLocalVector(x))));
47 }
48 
49 HypreParVector::HypreParVector(MPI_Comm comm, HYPRE_Int glob_size,
50  HYPRE_Int *col) : Vector()
51 {
52  x = hypre_ParVectorCreate(comm,glob_size,col);
53  hypre_ParVectorInitialize(x);
54  hypre_ParVectorSetPartitioningOwner(x,0);
55  // The data will be destroyed by hypre (this is the default)
56  hypre_ParVectorSetDataOwner(x,1);
57  hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
58  _SetDataAndSize_();
59  own_ParVector = 1;
60 }
61 
62 HypreParVector::HypreParVector(MPI_Comm comm, HYPRE_Int glob_size,
63  double *_data, HYPRE_Int *col) : Vector()
64 {
65  x = hypre_ParVectorCreate(comm,glob_size,col);
66  hypre_ParVectorSetDataOwner(x,1); // owns the seq vector
67  hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),0);
68  hypre_ParVectorSetPartitioningOwner(x,0);
69  double tmp = 0.0;
70  hypre_VectorData(hypre_ParVectorLocalVector(x)) = &tmp;
71  // If hypre_ParVectorLocalVector(x) and &tmp are non-NULL,
72  // hypre_ParVectorInitialize(x) does not allocate memory!
73  hypre_ParVectorInitialize(x);
74  // Set the internal data array to the one passed in
75  hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
76  _SetDataAndSize_();
77  own_ParVector = 1;
78 }
79 
81 {
82  x = hypre_ParVectorCreate(y.x -> comm, y.x -> global_size,
83  y.x -> partitioning);
84  hypre_ParVectorInitialize(x);
85  hypre_ParVectorSetPartitioningOwner(x,0);
86  hypre_ParVectorSetDataOwner(x,1);
87  hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
88  _SetDataAndSize_();
89  own_ParVector = 1;
90 }
91 
93  int transpose) : Vector()
94 {
95  if (!transpose)
96  {
97  x = hypre_ParVectorInDomainOf(const_cast<HypreParMatrix&>(A));
98  }
99  else
100  {
101  x = hypre_ParVectorInRangeOf(const_cast<HypreParMatrix&>(A));
102  }
103  _SetDataAndSize_();
104  own_ParVector = 1;
105 }
106 
108 {
109  x = (hypre_ParVector *) y;
110  _SetDataAndSize_();
111  own_ParVector = 0;
112 }
113 
115 {
116  x = hypre_ParVectorCreate(pfes->GetComm(), pfes->GlobalTrueVSize(),
117  pfes->GetTrueDofOffsets());
118  hypre_ParVectorInitialize(x);
119  hypre_ParVectorSetPartitioningOwner(x,0);
120  // The data will be destroyed by hypre (this is the default)
121  hypre_ParVectorSetDataOwner(x,1);
122  hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
123  _SetDataAndSize_();
124  own_ParVector = 1;
125 }
126 
128 {
129  hypre_Vector *hv = hypre_ParVectorToVectorAll(*this);
130  Vector *v = new Vector(hv->data, internal::to_int(hv->size));
131  v->MakeDataOwner();
132  hypre_SeqVectorSetDataOwner(hv,0);
133  hypre_SeqVectorDestroy(hv);
134  return v;
135 }
136 
138 {
139  hypre_ParVectorSetConstantValues(x,d);
140  return *this;
141 }
142 
144 {
145 #ifdef MFEM_DEBUG
146  if (size != y.Size())
147  {
148  mfem_error("HypreParVector::operator=");
149  }
150 #endif
151 
152  for (int i = 0; i < size; i++)
153  {
154  data[i] = y.data[i];
155  }
156  return *this;
157 }
158 
159 void HypreParVector::SetData(double *_data)
160 {
161  Vector::data = hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
162 }
163 
164 HYPRE_Int HypreParVector::Randomize(HYPRE_Int seed)
165 {
166  return hypre_ParVectorSetRandomValues(x,seed);
167 }
168 
169 void HypreParVector::Print(const char *fname) const
170 {
171  hypre_ParVectorPrint(x,fname);
172 }
173 
175 {
176  if (own_ParVector)
177  {
178  hypre_ParVectorDestroy(x);
179  }
180 }
181 
182 #ifdef MFEM_USE_SUNDIALS
183 
184 #ifndef SUNFALSE
185 #define SUNFALSE FALSE
186 #endif
187 
188 void HypreParVector::ToNVector(N_Vector &nv)
189 {
190  MFEM_ASSERT(nv && N_VGetVectorID(nv) == SUNDIALS_NVEC_PARHYP,
191  "invalid N_Vector");
192  N_VectorContent_ParHyp nv_c = (N_VectorContent_ParHyp)(nv->content);
193  MFEM_ASSERT(nv_c->own_parvector == SUNFALSE, "invalid N_Vector");
194  nv_c->local_length = x->local_vector->size;
195  nv_c->global_length = x->global_size;
196  nv_c->comm = x->comm;
197  nv_c->x = x;
198 }
199 
200 #endif // MFEM_USE_SUNDIALS
201 
202 
204 {
205  return hypre_ParVectorInnerProd(*x, *y);
206 }
207 
209 {
210  return hypre_ParVectorInnerProd(x, y);
211 }
212 
213 
214 double ParNormlp(const Vector &vec, double p, MPI_Comm comm)
215 {
216  double norm = 0.0;
217  if (p == 1.0)
218  {
219  double loc_norm = vec.Norml1();
220  MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
221  }
222  if (p == 2.0)
223  {
224  double loc_norm = vec*vec;
225  MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
226  norm = sqrt(norm);
227  }
228  if (p < std::numeric_limits<double>::infinity())
229  {
230  double sum = 0.0;
231  for (int i = 0; i < vec.Size(); i++)
232  {
233  sum += pow(fabs(vec(i)), p);
234  }
235  MPI_Allreduce(&sum, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
236  norm = pow(norm, 1.0/p);
237  }
238  else
239  {
240  double loc_norm = vec.Normlinf();
241  MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_MAX, comm);
242  }
243  return norm;
244 }
245 
246 
247 void HypreParMatrix::Init()
248 {
249  A = NULL;
250  X = Y = NULL;
251  diagOwner = offdOwner = colMapOwner = -1;
252  ParCSROwner = 1;
253 }
254 
256 {
257  Init();
258  height = width = 0;
259 }
260 
261 char HypreParMatrix::CopyCSR(SparseMatrix *csr, hypre_CSRMatrix *hypre_csr)
262 {
263  hypre_CSRMatrixData(hypre_csr) = csr->GetData();
264 #ifndef HYPRE_BIGINT
265  hypre_CSRMatrixI(hypre_csr) = csr->GetI();
266  hypre_CSRMatrixJ(hypre_csr) = csr->GetJ();
267  // Prevent hypre from destroying hypre_csr->{i,j,data}
268  return 0;
269 #else
270  hypre_CSRMatrixI(hypre_csr) =
271  DuplicateAs<HYPRE_Int>(csr->GetI(), csr->Height()+1);
272  hypre_CSRMatrixJ(hypre_csr) =
273  DuplicateAs<HYPRE_Int>(csr->GetJ(), csr->NumNonZeroElems());
274  // Prevent hypre from destroying hypre_csr->{i,j,data}, own {i,j}
275  return 1;
276 #endif
277 }
278 
279 char HypreParMatrix::CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr)
280 {
281  int nnz = bool_csr->Size_of_connections();
282  double *data = new double[nnz];
283  for (int i = 0; i < nnz; i++)
284  {
285  data[i] = 1.0;
286  }
287  hypre_CSRMatrixData(hypre_csr) = data;
288 #ifndef HYPRE_BIGINT
289  hypre_CSRMatrixI(hypre_csr) = bool_csr->GetI();
290  hypre_CSRMatrixJ(hypre_csr) = bool_csr->GetJ();
291  // Prevent hypre from destroying hypre_csr->{i,j,data}, own {data}
292  return 2;
293 #else
294  hypre_CSRMatrixI(hypre_csr) =
295  DuplicateAs<HYPRE_Int>(bool_csr->GetI(), bool_csr->Size()+1);
296  hypre_CSRMatrixJ(hypre_csr) =
297  DuplicateAs<HYPRE_Int>(bool_csr->GetJ(), nnz);
298  // Prevent hypre from destroying hypre_csr->{i,j,data}, own {i,j,data}
299  return 3;
300 #endif
301 }
302 
303 void HypreParMatrix::CopyCSR_J(hypre_CSRMatrix *hypre_csr, int *J)
304 {
305  HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(hypre_csr);
306  for (HYPRE_Int j = 0; j < nnz; j++)
307  {
308  J[j] = int(hypre_CSRMatrixJ(hypre_csr)[j]);
309  }
310 }
311 
312 // Square block-diagonal constructor
313 HypreParMatrix::HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size,
314  HYPRE_Int *row_starts, SparseMatrix *diag)
315  : Operator(diag->Height(), diag->Width())
316 {
317  Init();
318  A = hypre_ParCSRMatrixCreate(comm, glob_size, glob_size, row_starts,
319  row_starts, 0, diag->NumNonZeroElems(), 0);
320  hypre_ParCSRMatrixSetDataOwner(A,1);
321  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
322  hypre_ParCSRMatrixSetColStartsOwner(A,0);
323 
324  hypre_CSRMatrixSetDataOwner(A->diag,0);
325  diagOwner = CopyCSR(diag, A->diag);
326  hypre_CSRMatrixSetRownnz(A->diag);
327 
328  hypre_CSRMatrixSetDataOwner(A->offd,1);
329  hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->Height()+1);
330 
331  /* Don't need to call these, since they allocate memory only
332  if it was not already allocated */
333  // hypre_CSRMatrixInitialize(A->diag);
334  // hypre_ParCSRMatrixInitialize(A);
335 
336  hypre_ParCSRMatrixSetNumNonzeros(A);
337 
338  /* Make sure that the first entry in each row is the diagonal one. */
339  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
340 #ifdef HYPRE_BIGINT
341  CopyCSR_J(A->diag, diag->GetJ());
342 #endif
343 
344  hypre_MatvecCommPkgCreate(A);
345 }
346 
347 // Rectangular block-diagonal constructor
349  HYPRE_Int global_num_rows,
350  HYPRE_Int global_num_cols,
351  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
352  SparseMatrix *diag)
353  : Operator(diag->Height(), diag->Width())
354 {
355  Init();
356  A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
357  row_starts, col_starts,
358  0, diag->NumNonZeroElems(), 0);
359  hypre_ParCSRMatrixSetDataOwner(A,1);
360  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
361  hypre_ParCSRMatrixSetColStartsOwner(A,0);
362 
363  hypre_CSRMatrixSetDataOwner(A->diag,0);
364  diagOwner = CopyCSR(diag, A->diag);
365  hypre_CSRMatrixSetRownnz(A->diag);
366 
367  hypre_CSRMatrixSetDataOwner(A->offd,1);
368  hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->Height()+1);
369 
370  hypre_ParCSRMatrixSetNumNonzeros(A);
371 
372  /* Make sure that the first entry in each row is the diagonal one. */
373  if (row_starts == col_starts)
374  {
375  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
376 #ifdef HYPRE_BIGINT
377  CopyCSR_J(A->diag, diag->GetJ());
378 #endif
379  }
380 
381  hypre_MatvecCommPkgCreate(A);
382 }
383 
384 // General rectangular constructor with diagonal and off-diagonal
386  HYPRE_Int global_num_rows,
387  HYPRE_Int global_num_cols,
388  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
389  SparseMatrix *diag, SparseMatrix *offd,
390  HYPRE_Int *cmap)
391  : Operator(diag->Height(), diag->Width())
392 {
393  Init();
394  A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
395  row_starts, col_starts,
396  offd->Width(), diag->NumNonZeroElems(),
397  offd->NumNonZeroElems());
398  hypre_ParCSRMatrixSetDataOwner(A,1);
399  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
400  hypre_ParCSRMatrixSetColStartsOwner(A,0);
401 
402  hypre_CSRMatrixSetDataOwner(A->diag,0);
403  diagOwner = CopyCSR(diag, A->diag);
404  hypre_CSRMatrixSetRownnz(A->diag);
405 
406  hypre_CSRMatrixSetDataOwner(A->offd,0);
407  offdOwner = CopyCSR(offd, A->offd);
408  hypre_CSRMatrixSetRownnz(A->offd);
409 
410  hypre_ParCSRMatrixColMapOffd(A) = cmap;
411  // Prevent hypre from destroying A->col_map_offd
412  colMapOwner = 0;
413 
414  hypre_ParCSRMatrixSetNumNonzeros(A);
415 
416  /* Make sure that the first entry in each row is the diagonal one. */
417  if (row_starts == col_starts)
418  {
419  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
420 #ifdef HYPRE_BIGINT
421  CopyCSR_J(A->diag, diag->GetJ());
422 #endif
423  }
424 
425  hypre_MatvecCommPkgCreate(A);
426 }
427 
428 // General rectangular constructor with diagonal and off-diagonal
430  MPI_Comm comm,
431  HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
432  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
433  HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data,
434  HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data,
435  HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map)
436 {
437  Init();
438  A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
439  row_starts, col_starts, offd_num_cols, 0, 0);
440  hypre_ParCSRMatrixSetDataOwner(A,1);
441  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
442  hypre_ParCSRMatrixSetColStartsOwner(A,0);
443 
444  HYPRE_Int local_num_rows = hypre_CSRMatrixNumRows(A->diag);
445 
446  hypre_CSRMatrixSetDataOwner(A->diag,0);
447  hypre_CSRMatrixI(A->diag) = diag_i;
448  hypre_CSRMatrixJ(A->diag) = diag_j;
449  hypre_CSRMatrixData(A->diag) = diag_data;
450  hypre_CSRMatrixNumNonzeros(A->diag) = diag_i[local_num_rows];
451  hypre_CSRMatrixSetRownnz(A->diag);
452  // Prevent hypre from destroying A->diag->{i,j,data}, own A->diag->{i,j,data}
453  diagOwner = 3;
454 
455  hypre_CSRMatrixSetDataOwner(A->offd,0);
456  hypre_CSRMatrixI(A->offd) = offd_i;
457  hypre_CSRMatrixJ(A->offd) = offd_j;
458  hypre_CSRMatrixData(A->offd) = offd_data;
459  hypre_CSRMatrixNumNonzeros(A->offd) = offd_i[local_num_rows];
460  hypre_CSRMatrixSetRownnz(A->offd);
461  // Prevent hypre from destroying A->offd->{i,j,data}, own A->offd->{i,j,data}
462  offdOwner = 3;
463 
464  hypre_ParCSRMatrixColMapOffd(A) = offd_col_map;
465  // Prevent hypre from destroying A->col_map_offd, own A->col_map_offd
466  colMapOwner = 1;
467 
468  hypre_ParCSRMatrixSetNumNonzeros(A);
469 
470  /* Make sure that the first entry in each row is the diagonal one. */
471  if (row_starts == col_starts)
472  {
473  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
474  }
475 
476  hypre_MatvecCommPkgCreate(A);
477 
478  height = GetNumRows();
479  width = GetNumCols();
480 }
481 
482 // Constructor from a CSR matrix on rank 0
484  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
485  SparseMatrix *sm_a)
486 {
487  MFEM_ASSERT(sm_a != NULL, "invalid input");
488  MFEM_VERIFY(!HYPRE_AssumedPartitionCheck(),
489  "this method can not be used with assumed partition");
490 
491  Init();
492 
493  hypre_CSRMatrix *csr_a;
494  csr_a = hypre_CSRMatrixCreate(sm_a -> Height(), sm_a -> Width(),
495  sm_a -> NumNonZeroElems());
496 
497  hypre_CSRMatrixSetDataOwner(csr_a,0);
498  CopyCSR(sm_a, csr_a);
499  hypre_CSRMatrixSetRownnz(csr_a);
500 
501  A = hypre_CSRMatrixToParCSRMatrix(comm, csr_a, row_starts, col_starts);
502 
503 #ifdef HYPRE_BIGINT
504  delete [] hypre_CSRMatrixI(csr_a);
505  delete [] hypre_CSRMatrixJ(csr_a);
506 #endif
507  hypre_CSRMatrixI(csr_a) = NULL;
508  hypre_CSRMatrixDestroy(csr_a);
509 
510  height = GetNumRows();
511  width = GetNumCols();
512 
513  /* Make sure that the first entry in each row is the diagonal one. */
514  if (row_starts == col_starts)
515  {
516  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
517  }
518 
519  hypre_MatvecCommPkgCreate(A);
520 }
521 
522 // Boolean, rectangular, block-diagonal constructor
524  HYPRE_Int global_num_rows,
525  HYPRE_Int global_num_cols,
526  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
527  Table *diag)
528 {
529  Init();
530  int nnz = diag->Size_of_connections();
531  A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
532  row_starts, col_starts, 0, nnz, 0);
533  hypre_ParCSRMatrixSetDataOwner(A,1);
534  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
535  hypre_ParCSRMatrixSetColStartsOwner(A,0);
536 
537  hypre_CSRMatrixSetDataOwner(A->diag,0);
538  diagOwner = CopyBoolCSR(diag, A->diag);
539  hypre_CSRMatrixSetRownnz(A->diag);
540 
541  hypre_CSRMatrixSetDataOwner(A->offd,1);
542  hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->Size()+1);
543 
544  hypre_ParCSRMatrixSetNumNonzeros(A);
545 
546  /* Make sure that the first entry in each row is the diagonal one. */
547  if (row_starts == col_starts)
548  {
549  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
550 #ifdef HYPRE_BIGINT
551  CopyCSR_J(A->diag, diag->GetJ());
552 #endif
553  }
554 
555  hypre_MatvecCommPkgCreate(A);
556 
557  height = GetNumRows();
558  width = GetNumCols();
559 }
560 
561 // Boolean, general rectangular constructor with diagonal and off-diagonal
562 HypreParMatrix::HypreParMatrix(MPI_Comm comm, int id, int np,
563  HYPRE_Int *row, HYPRE_Int *col,
564  HYPRE_Int *i_diag, HYPRE_Int *j_diag,
565  HYPRE_Int *i_offd, HYPRE_Int *j_offd,
566  HYPRE_Int *cmap, HYPRE_Int cmap_size)
567 {
568  HYPRE_Int diag_nnz, offd_nnz;
569 
570  Init();
571  if (HYPRE_AssumedPartitionCheck())
572  {
573  diag_nnz = i_diag[row[1]-row[0]];
574  offd_nnz = i_offd[row[1]-row[0]];
575 
576  A = hypre_ParCSRMatrixCreate(comm, row[2], col[2], row, col,
577  cmap_size, diag_nnz, offd_nnz);
578  }
579  else
580  {
581  diag_nnz = i_diag[row[id+1]-row[id]];
582  offd_nnz = i_offd[row[id+1]-row[id]];
583 
584  A = hypre_ParCSRMatrixCreate(comm, row[np], col[np], row, col,
585  cmap_size, diag_nnz, offd_nnz);
586  }
587 
588  hypre_ParCSRMatrixSetDataOwner(A,1);
589  hypre_ParCSRMatrixSetRowStartsOwner(A,0);
590  hypre_ParCSRMatrixSetColStartsOwner(A,0);
591 
592  HYPRE_Int i;
593 
594  double *a_diag = new double[diag_nnz];
595  for (i = 0; i < diag_nnz; i++)
596  {
597  a_diag[i] = 1.0;
598  }
599 
600  double *a_offd = new double[offd_nnz];
601  for (i = 0; i < offd_nnz; i++)
602  {
603  a_offd[i] = 1.0;
604  }
605 
606  hypre_CSRMatrixSetDataOwner(A->diag,0);
607  hypre_CSRMatrixI(A->diag) = i_diag;
608  hypre_CSRMatrixJ(A->diag) = j_diag;
609  hypre_CSRMatrixData(A->diag) = a_diag;
610  hypre_CSRMatrixSetRownnz(A->diag);
611  // Prevent hypre from destroying A->diag->{i,j,data}, own A->diag->{i,j,data}
612  diagOwner = 3;
613 
614  hypre_CSRMatrixSetDataOwner(A->offd,0);
615  hypre_CSRMatrixI(A->offd) = i_offd;
616  hypre_CSRMatrixJ(A->offd) = j_offd;
617  hypre_CSRMatrixData(A->offd) = a_offd;
618  hypre_CSRMatrixSetRownnz(A->offd);
619  // Prevent hypre from destroying A->offd->{i,j,data}, own A->offd->{i,j,data}
620  offdOwner = 3;
621 
622  hypre_ParCSRMatrixColMapOffd(A) = cmap;
623  // Prevent hypre from destroying A->col_map_offd, own A->col_map_offd
624  colMapOwner = 1;
625 
626  hypre_ParCSRMatrixSetNumNonzeros(A);
627 
628  /* Make sure that the first entry in each row is the diagonal one. */
629  if (row == col)
630  {
631  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
632  }
633 
634  hypre_MatvecCommPkgCreate(A);
635 
636  height = GetNumRows();
637  width = GetNumCols();
638 }
639 
640 // General rectangular constructor with diagonal and off-diagonal constructed
641 // from a CSR matrix that contains both diagonal and off-diagonal blocks
642 HypreParMatrix::HypreParMatrix(MPI_Comm comm, int nrows, HYPRE_Int glob_nrows,
643  HYPRE_Int glob_ncols, int *I, HYPRE_Int *J,
644  double *data, HYPRE_Int *rows, HYPRE_Int *cols)
645 {
646  Init();
647 
648  // Determine partitioning size, and my column start and end
649  int part_size;
650  HYPRE_Int my_col_start, my_col_end; // my range: [my_col_start, my_col_end)
651  if (HYPRE_AssumedPartitionCheck())
652  {
653  part_size = 2;
654  my_col_start = cols[0];
655  my_col_end = cols[1];
656  }
657  else
658  {
659  int myid;
660  MPI_Comm_rank(comm, &myid);
661  MPI_Comm_size(comm, &part_size);
662  part_size++;
663  my_col_start = cols[myid];
664  my_col_end = cols[myid+1];
665  }
666 
667  // Copy in the row and column partitionings
668  HYPRE_Int *row_starts, *col_starts;
669  if (rows == cols)
670  {
671  row_starts = col_starts = hypre_TAlloc(HYPRE_Int, part_size);
672  for (int i = 0; i < part_size; i++)
673  {
674  row_starts[i] = rows[i];
675  }
676  }
677  else
678  {
679  row_starts = hypre_TAlloc(HYPRE_Int, part_size);
680  col_starts = hypre_TAlloc(HYPRE_Int, part_size);
681  for (int i = 0; i < part_size; i++)
682  {
683  row_starts[i] = rows[i];
684  col_starts[i] = cols[i];
685  }
686  }
687 
688  // Create a map for the off-diagonal indices - global to local. Count the
689  // number of diagonal and off-diagonal entries.
690  HYPRE_Int diag_nnz = 0, offd_nnz = 0, offd_num_cols = 0;
691  map<HYPRE_Int, HYPRE_Int> offd_map;
692  for (HYPRE_Int j = 0, loc_nnz = I[nrows]; j < loc_nnz; j++)
693  {
694  HYPRE_Int glob_col = J[j];
695  if (my_col_start <= glob_col && glob_col < my_col_end)
696  {
697  diag_nnz++;
698  }
699  else
700  {
701  offd_map.insert(pair<const HYPRE_Int, HYPRE_Int>(glob_col, -1));
702  offd_nnz++;
703  }
704  }
705  // count the number of columns in the off-diagonal and set the local indices
706  for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
707  it != offd_map.end(); ++it)
708  {
709  it->second = offd_num_cols++;
710  }
711 
712  // construct the global ParCSR matrix
713  A = hypre_ParCSRMatrixCreate(comm, glob_nrows, glob_ncols,
714  row_starts, col_starts, offd_num_cols,
715  diag_nnz, offd_nnz);
716  hypre_ParCSRMatrixInitialize(A);
717 
718  HYPRE_Int *diag_i, *diag_j, *offd_i, *offd_j, *offd_col_map;
719  double *diag_data, *offd_data;
720  diag_i = A->diag->i;
721  diag_j = A->diag->j;
722  diag_data = A->diag->data;
723  offd_i = A->offd->i;
724  offd_j = A->offd->j;
725  offd_data = A->offd->data;
726  offd_col_map = A->col_map_offd;
727 
728  diag_nnz = offd_nnz = 0;
729  for (HYPRE_Int i = 0, j = 0; i < nrows; i++)
730  {
731  diag_i[i] = diag_nnz;
732  offd_i[i] = offd_nnz;
733  for (HYPRE_Int j_end = I[i+1]; j < j_end; j++)
734  {
735  HYPRE_Int glob_col = J[j];
736  if (my_col_start <= glob_col && glob_col < my_col_end)
737  {
738  diag_j[diag_nnz] = glob_col - my_col_start;
739  diag_data[diag_nnz] = data[j];
740  diag_nnz++;
741  }
742  else
743  {
744  offd_j[offd_nnz] = offd_map[glob_col];
745  offd_data[offd_nnz] = data[j];
746  offd_nnz++;
747  }
748  }
749  }
750  diag_i[nrows] = diag_nnz;
751  offd_i[nrows] = offd_nnz;
752  for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
753  it != offd_map.end(); ++it)
754  {
755  offd_col_map[it->second] = it->first;
756  }
757 
758  hypre_ParCSRMatrixSetNumNonzeros(A);
759  /* Make sure that the first entry in each row is the diagonal one. */
760  if (row_starts == col_starts)
761  {
762  hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
763  }
764  hypre_MatvecCommPkgCreate(A);
765 
766  height = GetNumRows();
767  width = GetNumCols();
768 }
769 
771 {
772  Destroy();
773  Init();
774  A = master.A;
775  ParCSROwner = 0;
776  height = master.GetNumRows();
777  width = master.GetNumCols();
778 }
779 
780 hypre_ParCSRMatrix* HypreParMatrix::StealData()
781 {
782  // Only safe when (diagOwner == -1 && offdOwner == -1 && colMapOwner == -1)
783  // Otherwise, there may be memory leaks or hypre may destroy arrays allocated
784  // with operator new.
785  MFEM_ASSERT(diagOwner == -1 && offdOwner == -1 && colMapOwner == -1, "");
786  MFEM_ASSERT(ParCSROwner, "");
787  hypre_ParCSRMatrix *R = A;
788  A = NULL;
789  Destroy();
790  Init();
791  return R;
792 }
793 
795 {
796  if (!A || hypre_ParCSRMatrixOwnsRowStarts(A) ||
797  (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
798  hypre_ParCSRMatrixOwnsColStarts(A)))
799  {
800  return;
801  }
802 
803  int row_starts_size;
804  if (HYPRE_AssumedPartitionCheck())
805  {
806  row_starts_size = 2;
807  }
808  else
809  {
810  MPI_Comm_size(hypre_ParCSRMatrixComm(A), &row_starts_size);
811  row_starts_size++; // num_proc + 1
812  }
813 
814  HYPRE_Int *old_row_starts = hypre_ParCSRMatrixRowStarts(A);
815  HYPRE_Int *new_row_starts = hypre_CTAlloc(HYPRE_Int, row_starts_size);
816  for (int i = 0; i < row_starts_size; i++)
817  {
818  new_row_starts[i] = old_row_starts[i];
819  }
820 
821  hypre_ParCSRMatrixRowStarts(A) = new_row_starts;
822  hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
823 
824  if (hypre_ParCSRMatrixColStarts(A) == old_row_starts)
825  {
826  hypre_ParCSRMatrixColStarts(A) = new_row_starts;
827  hypre_ParCSRMatrixOwnsColStarts(A) = 0;
828  }
829 }
830 
832 {
833  if (!A || hypre_ParCSRMatrixOwnsColStarts(A) ||
834  (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
835  hypre_ParCSRMatrixOwnsRowStarts(A)))
836  {
837  return;
838  }
839 
840  int col_starts_size;
841  if (HYPRE_AssumedPartitionCheck())
842  {
843  col_starts_size = 2;
844  }
845  else
846  {
847  MPI_Comm_size(hypre_ParCSRMatrixComm(A), &col_starts_size);
848  col_starts_size++; // num_proc + 1
849  }
850 
851  HYPRE_Int *old_col_starts = hypre_ParCSRMatrixColStarts(A);
852  HYPRE_Int *new_col_starts = hypre_CTAlloc(HYPRE_Int, col_starts_size);
853  for (int i = 0; i < col_starts_size; i++)
854  {
855  new_col_starts[i] = old_col_starts[i];
856  }
857 
858  hypre_ParCSRMatrixColStarts(A) = new_col_starts;
859 
860  if (hypre_ParCSRMatrixRowStarts(A) == old_col_starts)
861  {
862  hypre_ParCSRMatrixRowStarts(A) = new_col_starts;
863  hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
864  hypre_ParCSRMatrixOwnsColStarts(A) = 0;
865  }
866  else
867  {
868  hypre_ParCSRMatrixOwnsColStarts(A) = 1;
869  }
870 }
871 
873 {
874  int size = Height();
875  diag.SetSize(size);
876  for (int j = 0; j < size; j++)
877  {
878  diag(j) = A->diag->data[A->diag->i[j]];
879  MFEM_ASSERT(A->diag->j[A->diag->i[j]] == j,
880  "the first entry in each row must be the diagonal one");
881  }
882 }
883 
884 static void MakeWrapper(const hypre_CSRMatrix *mat, SparseMatrix &wrapper)
885 {
886  HYPRE_Int nr = hypre_CSRMatrixNumRows(mat);
887  HYPRE_Int nc = hypre_CSRMatrixNumCols(mat);
888 #ifndef HYPRE_BIGINT
889  SparseMatrix tmp(hypre_CSRMatrixI(mat),
890  hypre_CSRMatrixJ(mat),
891  hypre_CSRMatrixData(mat),
892  nr, nc, false, false, false);
893 #else
894  HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(mat);
895  SparseMatrix tmp(DuplicateAs<int>(hypre_CSRMatrixI(mat), nr+1),
896  DuplicateAs<int>(hypre_CSRMatrixJ(mat), nnz),
897  hypre_CSRMatrixData(mat),
898  nr, nc, true, false, false);
899 #endif
900  wrapper.Swap(tmp);
901 }
902 
904 {
905  MakeWrapper(A->diag, diag);
906 }
907 
908 void HypreParMatrix::GetOffd(SparseMatrix &offd, HYPRE_Int* &cmap) const
909 {
910  MakeWrapper(A->offd, offd);
911  cmap = A->col_map_offd;
912 }
913 
915  bool interleaved_rows,
916  bool interleaved_cols) const
917 {
918  int nr = blocks.NumRows();
919  int nc = blocks.NumCols();
920 
921  hypre_ParCSRMatrix **hypre_blocks = new hypre_ParCSRMatrix*[nr * nc];
922  internal::hypre_ParCSRMatrixSplit(A, nr, nc, hypre_blocks,
923  interleaved_rows, interleaved_cols);
924 
925  for (int i = 0; i < nr; i++)
926  {
927  for (int j = 0; j < nc; j++)
928  {
929  blocks[i][j] = new HypreParMatrix(hypre_blocks[i*nc + j]);
930  }
931  }
932 
933  delete [] hypre_blocks;
934 }
935 
937 {
938  hypre_ParCSRMatrix * At;
939  hypre_ParCSRMatrixTranspose(A, &At, 1);
940  hypre_ParCSRMatrixSetNumNonzeros(At);
941 
942  hypre_MatvecCommPkgCreate(At);
943 
944  return new HypreParMatrix(At);
945 }
946 
948  double a, double b)
949 {
950  return hypre_ParCSRMatrixMatvec(a, A, x, b, y);
951 }
952 
953 void HypreParMatrix::Mult(double a, const Vector &x, double b, Vector &y) const
954 {
955  MFEM_ASSERT(x.Size() == Width(), "invalid x.Size() = " << x.Size()
956  << ", expected size = " << Width());
957  MFEM_ASSERT(y.Size() == Height(), "invalid y.Size() = " << y.Size()
958  << ", expected size = " << Height());
959 
960  if (X == NULL)
961  {
962  X = new HypreParVector(A->comm,
964  x.GetData(),
965  GetColStarts());
966  Y = new HypreParVector(A->comm,
968  y.GetData(),
969  GetRowStarts());
970  }
971  else
972  {
973  X->SetData(x.GetData());
974  Y->SetData(y.GetData());
975  }
976 
977  hypre_ParCSRMatrixMatvec(a, A, *X, b, *Y);
978 }
979 
980 void HypreParMatrix::MultTranspose(double a, const Vector &x,
981  double b, Vector &y) const
982 {
983  MFEM_ASSERT(x.Size() == Height(), "invalid x.Size() = " << x.Size()
984  << ", expected size = " << Height());
985  MFEM_ASSERT(y.Size() == Width(), "invalid y.Size() = " << y.Size()
986  << ", expected size = " << Width());
987 
988  // Note: x has the dimensions of Y (height), and
989  // y has the dimensions of X (width)
990  if (X == NULL)
991  {
992  X = new HypreParVector(A->comm,
994  y.GetData(),
995  GetColStarts());
996  Y = new HypreParVector(A->comm,
998  x.GetData(),
999  GetRowStarts());
1000  }
1001  else
1002  {
1003  X->SetData(y.GetData());
1004  Y->SetData(x.GetData());
1005  }
1006 
1007  hypre_ParCSRMatrixMatvecT(a, A, *Y, b, *X);
1008 }
1009 
1010 HYPRE_Int HypreParMatrix::Mult(HYPRE_ParVector x, HYPRE_ParVector y,
1011  double a, double b)
1012 {
1013  return hypre_ParCSRMatrixMatvec(a, A, (hypre_ParVector *) x, b,
1014  (hypre_ParVector *) y);
1015 }
1016 
1018  double a, double b)
1019 {
1020  return hypre_ParCSRMatrixMatvecT(a, A, x, b, y);
1021 }
1022 
1024  HYPRE_Int* row_starts) const
1025 {
1026  const bool assumed_partition = HYPRE_AssumedPartitionCheck();
1027  if (row_starts == NULL)
1028  {
1029  row_starts = hypre_ParCSRMatrixRowStarts(A);
1030  MFEM_VERIFY(D.Height() == hypre_CSRMatrixNumRows(A->diag),
1031  "the matrix D is NOT compatible with the row starts of"
1032  " this HypreParMatrix, row_starts must be given.");
1033  }
1034  else
1035  {
1036  int offset;
1037  if (assumed_partition)
1038  {
1039  offset = 0;
1040  }
1041  else
1042  {
1043  MPI_Comm_rank(GetComm(), &offset);
1044  }
1045  int local_num_rows = row_starts[offset+1]-row_starts[offset];
1046  MFEM_VERIFY(local_num_rows == D.Height(), "the number of rows in D is "
1047  " not compatible with the given row_starts");
1048  }
1049  // D.Width() will be checked for compatibility by the SparseMatrix
1050  // multiplication function, mfem::Mult(), called below.
1051 
1052  int part_size;
1053  HYPRE_Int global_num_rows;
1054  if (assumed_partition)
1055  {
1056  part_size = 2;
1057  global_num_rows = row_starts[2];
1058  // Here, we use row_starts[2], so row_starts must come from the methods
1059  // GetDofOffsets/GetTrueDofOffsets of ParFiniteElementSpace (HYPRE's
1060  // partitions have only 2 entries).
1061  }
1062  else
1063  {
1064  MPI_Comm_size(GetComm(), &part_size);
1065  global_num_rows = row_starts[part_size];
1066  part_size++;
1067  }
1068 
1069  HYPRE_Int *col_starts = hypre_ParCSRMatrixColStarts(A);
1070  HYPRE_Int *col_map_offd;
1071 
1072  // get the diag and offd blocks as SparseMatrix wrappers
1073  SparseMatrix A_diag, A_offd;
1074  GetDiag(A_diag);
1075  GetOffd(A_offd, col_map_offd);
1076 
1077  // multiply the blocks with D and create a new HypreParMatrix
1078  SparseMatrix* DA_diag = mfem::Mult(D, A_diag);
1079  SparseMatrix* DA_offd = mfem::Mult(D, A_offd);
1080 
1081  HypreParMatrix* DA =
1082  new HypreParMatrix(GetComm(),
1083  global_num_rows, hypre_ParCSRMatrixGlobalNumCols(A),
1084  DuplicateAs<HYPRE_Int>(row_starts, part_size, false),
1085  DuplicateAs<HYPRE_Int>(col_starts, part_size, false),
1086  DA_diag, DA_offd,
1087  DuplicateAs<HYPRE_Int>(col_map_offd, A_offd.Width()));
1088 
1089  // When HYPRE_BIGINT is defined, we want DA_{diag,offd} to delete their I and
1090  // J arrays but not their data arrays; when HYPRE_BIGINT is not defined, we
1091  // don't want DA_{diag,offd} to delete anything.
1092 #ifndef HYPRE_BIGINT
1093  DA_diag->LoseData();
1094  DA_offd->LoseData();
1095 #else
1096  DA_diag->SetDataOwner(false);
1097  DA_offd->SetDataOwner(false);
1098 #endif
1099 
1100  delete DA_diag;
1101  delete DA_offd;
1102 
1103  hypre_ParCSRMatrixSetRowStartsOwner(DA->A, 1);
1104  hypre_ParCSRMatrixSetColStartsOwner(DA->A, 1);
1105 
1106  DA->diagOwner = DA->offdOwner = 3;
1107  DA->colMapOwner = 1;
1108 
1109  return DA;
1110 }
1111 
1113 {
1114  if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1115  {
1116  mfem_error("Row does not match");
1117  }
1118 
1119  if (hypre_CSRMatrixNumRows(A->diag) != diag.Size())
1120  {
1121  mfem_error("Note the Vector diag is not of compatible dimensions with A\n");
1122  }
1123 
1124  int size = Height();
1125  double *Adiag_data = hypre_CSRMatrixData(A->diag);
1126  HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1127 
1128 
1129  double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1130  HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1131  double val;
1132  HYPRE_Int jj;
1133  for (int i(0); i < size; ++i)
1134  {
1135  val = diag[i];
1136  for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1137  {
1138  Adiag_data[jj] *= val;
1139  }
1140  for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1141  {
1142  Aoffd_data[jj] *= val;
1143  }
1144  }
1145 }
1146 
1148 {
1149  if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1150  {
1151  mfem_error("Row does not match");
1152  }
1153 
1154  if (hypre_CSRMatrixNumRows(A->diag) != diag.Size())
1155  {
1156  mfem_error("Note the Vector diag is not of compatible dimensions with A\n");
1157  }
1158 
1159  int size = Height();
1160  double *Adiag_data = hypre_CSRMatrixData(A->diag);
1161  HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1162 
1163 
1164  double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1165  HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1166  double val;
1167  HYPRE_Int jj;
1168  for (int i(0); i < size; ++i)
1169  {
1170 #ifdef MFEM_DEBUG
1171  if (0.0 == diag(i))
1172  {
1173  mfem_error("HypreParMatrix::InvDiagScale : Division by 0");
1174  }
1175 #endif
1176  val = 1./diag(i);
1177  for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1178  {
1179  Adiag_data[jj] *= val;
1180  }
1181  for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1182  {
1183  Aoffd_data[jj] *= val;
1184  }
1185  }
1186 }
1187 
1189 {
1190  if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1191  {
1192  mfem_error("Row does not match");
1193  }
1194 
1195  HYPRE_Int size=hypre_CSRMatrixNumRows(A->diag);
1196  HYPRE_Int jj;
1197 
1198  double *Adiag_data = hypre_CSRMatrixData(A->diag);
1199  HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1200  for (jj = 0; jj < Adiag_i[size]; ++jj)
1201  {
1202  Adiag_data[jj] *= s;
1203  }
1204 
1205  double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1206  HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1207  for (jj = 0; jj < Aoffd_i[size]; ++jj)
1208  {
1209  Aoffd_data[jj] *= s;
1210  }
1211 }
1212 
1213 static void get_sorted_rows_cols(const Array<int> &rows_cols,
1214  Array<HYPRE_Int> &hypre_sorted)
1215 {
1216  hypre_sorted.SetSize(rows_cols.Size());
1217  bool sorted = true;
1218  for (int i = 0; i < rows_cols.Size(); i++)
1219  {
1220  hypre_sorted[i] = rows_cols[i];
1221  if (i && rows_cols[i-1] > rows_cols[i]) { sorted = false; }
1222  }
1223  if (!sorted) { hypre_sorted.Sort(); }
1224 }
1225 
1226 void HypreParMatrix::Threshold(double threshold)
1227 {
1228  int ierr = 0;
1229 
1230  MPI_Comm comm;
1231  hypre_CSRMatrix * csr_A;
1232  hypre_CSRMatrix * csr_A_wo_z;
1233  hypre_ParCSRMatrix * parcsr_A_ptr;
1234  HYPRE_Int * row_starts = NULL; HYPRE_Int * col_starts = NULL;
1235  HYPRE_Int row_start = -1; HYPRE_Int row_end = -1;
1236  HYPRE_Int col_start = -1; HYPRE_Int col_end = -1;
1237 
1238  comm = hypre_ParCSRMatrixComm(A);
1239 
1240  ierr += hypre_ParCSRMatrixGetLocalRange(A,
1241  &row_start,&row_end,
1242  &col_start,&col_end );
1243 
1244  row_starts = hypre_ParCSRMatrixRowStarts(A);
1245  col_starts = hypre_ParCSRMatrixColStarts(A);
1246 
1247  bool old_owns_row = hypre_ParCSRMatrixOwnsRowStarts(A);
1248  bool old_owns_col = hypre_ParCSRMatrixOwnsColStarts(A);
1249  HYPRE_Int global_num_rows = hypre_ParCSRMatrixGlobalNumRows(A);
1250  HYPRE_Int global_num_cols = hypre_ParCSRMatrixGlobalNumCols(A);
1251  parcsr_A_ptr = hypre_ParCSRMatrixCreate(comm, global_num_rows,
1252  global_num_cols,
1253  row_starts, col_starts,
1254  0, 0, 0);
1255  hypre_ParCSRMatrixOwnsRowStarts(parcsr_A_ptr) = old_owns_row;
1256  hypre_ParCSRMatrixOwnsColStarts(parcsr_A_ptr) = old_owns_col;
1257 
1258  csr_A = hypre_MergeDiagAndOffd(A);
1259 
1260  csr_A_wo_z = hypre_CSRMatrixDeleteZeros(csr_A,threshold);
1261 
1262  /* hypre_CSRMatrixDeleteZeros will return a NULL pointer rather than a usable
1263  CSR matrix if it finds no non-zeros */
1264  if (csr_A_wo_z == NULL)
1265  {
1266  csr_A_wo_z = csr_A;
1267  }
1268  else
1269  {
1270  ierr += hypre_CSRMatrixDestroy(csr_A);
1271  }
1272 
1273  /* FIXME: GenerateDiagAndOffd() uses an int array of size equal to the number
1274  of columns in csr_A_wo_z which is the global number of columns in A. This
1275  does not scale well. */
1276  ierr += GenerateDiagAndOffd(csr_A_wo_z,parcsr_A_ptr,
1277  col_start,col_end);
1278 
1279  ierr += hypre_CSRMatrixDestroy(csr_A_wo_z);
1280 
1281  ierr += hypre_ParCSRMatrixDestroy(A);
1282 
1283  A = parcsr_A_ptr;
1284 }
1285 
1287  const HypreParVector &X,
1288  HypreParVector &B)
1289 {
1290  Array<HYPRE_Int> rc_sorted;
1291  get_sorted_rows_cols(rows_cols, rc_sorted);
1292 
1293  internal::hypre_ParCSRMatrixEliminateAXB(
1294  A, rc_sorted.Size(), rc_sorted.GetData(), X, B);
1295 }
1296 
1298 {
1299  Array<HYPRE_Int> rc_sorted;
1300  get_sorted_rows_cols(rows_cols, rc_sorted);
1301 
1302  hypre_ParCSRMatrix* Ae;
1303  internal::hypre_ParCSRMatrixEliminateAAe(
1304  A, &Ae, rc_sorted.Size(), rc_sorted.GetData());
1305 
1306  return new HypreParMatrix(Ae);
1307 }
1308 
1309 void HypreParMatrix::Print(const char *fname, HYPRE_Int offi, HYPRE_Int offj)
1310 {
1311  hypre_ParCSRMatrixPrintIJ(A,offi,offj,fname);
1312 }
1313 
1314 void HypreParMatrix::Read(MPI_Comm comm, const char *fname)
1315 {
1316  Destroy();
1317  Init();
1318 
1319  HYPRE_Int base_i, base_j;
1320  hypre_ParCSRMatrixReadIJ(comm, fname, &base_i, &base_j, &A);
1321  hypre_ParCSRMatrixSetNumNonzeros(A);
1322 
1323  hypre_MatvecCommPkgCreate(A);
1324 
1325  height = GetNumRows();
1326  width = GetNumCols();
1327 }
1328 
1329 void HypreParMatrix::Read_IJMatrix(MPI_Comm comm, const char *fname)
1330 {
1331  Destroy();
1332  Init();
1333 
1334  HYPRE_IJMatrix A_ij;
1335  HYPRE_IJMatrixRead(fname, comm, 5555, &A_ij); // HYPRE_PARCSR = 5555
1336 
1337  HYPRE_ParCSRMatrix A_parcsr;
1338  HYPRE_IJMatrixGetObject(A_ij, (void**) &A_parcsr);
1339 
1340  A = (hypre_ParCSRMatrix*)A_parcsr;
1341 
1342  hypre_ParCSRMatrixSetNumNonzeros(A);
1343 
1344  hypre_MatvecCommPkgCreate(A);
1345 
1346  height = GetNumRows();
1347  width = GetNumCols();
1348 }
1349 
1350 void HypreParMatrix::PrintCommPkg(std::ostream &out) const
1351 {
1352  hypre_ParCSRCommPkg *comm_pkg = A->comm_pkg;
1353  MPI_Comm comm = A->comm;
1354  char c = '\0';
1355  const int tag = 46801;
1356  int myid, nproc;
1357  MPI_Comm_rank(comm, &myid);
1358  MPI_Comm_size(comm, &nproc);
1359 
1360  if (myid != 0)
1361  {
1362  MPI_Recv(&c, 1, MPI_CHAR, myid-1, tag, comm, MPI_STATUS_IGNORE);
1363  }
1364  else
1365  {
1366  out << "\nHypreParMatrix: hypre_ParCSRCommPkg:\n";
1367  }
1368  out << "Rank " << myid << ":\n"
1369  " number of sends = " << comm_pkg->num_sends <<
1370  " (" << sizeof(double)*comm_pkg->send_map_starts[comm_pkg->num_sends] <<
1371  " bytes)\n"
1372  " number of recvs = " << comm_pkg->num_recvs <<
1373  " (" << sizeof(double)*comm_pkg->recv_vec_starts[comm_pkg->num_recvs] <<
1374  " bytes)\n";
1375  if (myid != nproc-1)
1376  {
1377  out << std::flush;
1378  MPI_Send(&c, 1, MPI_CHAR, myid+1, tag, comm);
1379  }
1380  else
1381  {
1382  out << std::endl;
1383  }
1384  MPI_Barrier(comm);
1385 }
1386 
1387 void HypreParMatrix::Destroy()
1388 {
1389  if ( X != NULL ) { delete X; }
1390  if ( Y != NULL ) { delete Y; }
1391 
1392  if (A == NULL) { return; }
1393 
1394  if (diagOwner >= 0)
1395  {
1396  if (diagOwner & 1)
1397  {
1398  delete [] hypre_CSRMatrixI(A->diag);
1399  delete [] hypre_CSRMatrixJ(A->diag);
1400  }
1401  hypre_CSRMatrixI(A->diag) = NULL;
1402  hypre_CSRMatrixJ(A->diag) = NULL;
1403  if (diagOwner & 2)
1404  {
1405  delete [] hypre_CSRMatrixData(A->diag);
1406  }
1407  hypre_CSRMatrixData(A->diag) = NULL;
1408  }
1409  if (offdOwner >= 0)
1410  {
1411  if (offdOwner & 1)
1412  {
1413  delete [] hypre_CSRMatrixI(A->offd);
1414  delete [] hypre_CSRMatrixJ(A->offd);
1415  }
1416  hypre_CSRMatrixI(A->offd) = NULL;
1417  hypre_CSRMatrixJ(A->offd) = NULL;
1418  if (offdOwner & 2)
1419  {
1420  delete [] hypre_CSRMatrixData(A->offd);
1421  }
1422  hypre_CSRMatrixData(A->offd) = NULL;
1423  }
1424  if (colMapOwner >= 0)
1425  {
1426  if (colMapOwner & 1)
1427  {
1428  delete [] hypre_ParCSRMatrixColMapOffd(A);
1429  }
1430  hypre_ParCSRMatrixColMapOffd(A) = NULL;
1431  }
1432 
1433  if (ParCSROwner)
1434  {
1435  hypre_ParCSRMatrixDestroy(A);
1436  }
1437 }
1438 
1440  double beta, const HypreParMatrix &B)
1441 {
1442  hypre_ParCSRMatrix *C_hypre =
1443  internal::hypre_ParCSRMatrixAdd(const_cast<HypreParMatrix &>(A),
1444  const_cast<HypreParMatrix &>(B));
1445  MFEM_VERIFY(C_hypre, "error in hypre_ParCSRMatrixAdd");
1446 
1447  hypre_MatvecCommPkgCreate(C_hypre);
1448  HypreParMatrix *C = new HypreParMatrix(C_hypre);
1449  *C = 0.0;
1450  C->Add(alpha, A);
1451  C->Add(beta, B);
1452 
1453  return C;
1454 }
1455 
1457 {
1458  hypre_ParCSRMatrix * ab;
1459  ab = hypre_ParMatmul(*A,*B);
1460  hypre_ParCSRMatrixSetNumNonzeros(ab);
1461 
1462  hypre_MatvecCommPkgCreate(ab);
1463 
1464  return new HypreParMatrix(ab);
1465 }
1466 
1468 {
1469  hypre_ParCSRMatrix * C = internal::hypre_ParCSRMatrixAdd(*A,*B);
1470 
1471  hypre_MatvecCommPkgCreate(C);
1472 
1473  return new HypreParMatrix(C);
1474 }
1475 
1477 {
1478  HYPRE_Int P_owns_its_col_starts =
1479  hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1480 
1481  hypre_ParCSRMatrix * rap;
1482  hypre_BoomerAMGBuildCoarseOperator(*P,*A,*P,&rap);
1483  hypre_ParCSRMatrixSetNumNonzeros(rap);
1484  // hypre_MatvecCommPkgCreate(rap);
1485 
1486  /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
1487  from P (even if it does not own them)! */
1488  hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1489  hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1490 
1491  if (P_owns_its_col_starts)
1492  {
1493  hypre_ParCSRMatrixSetColStartsOwner(*P, 1);
1494  }
1495 
1496  return new HypreParMatrix(rap);
1497 }
1498 
1500  const HypreParMatrix *P)
1501 {
1502  HYPRE_Int P_owns_its_col_starts =
1503  hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1504  HYPRE_Int Rt_owns_its_col_starts =
1505  hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*Rt));
1506 
1507  hypre_ParCSRMatrix * rap;
1508  hypre_BoomerAMGBuildCoarseOperator(*Rt,*A,*P,&rap);
1509 
1510  hypre_ParCSRMatrixSetNumNonzeros(rap);
1511  // hypre_MatvecCommPkgCreate(rap);
1512 
1513  /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
1514  from Rt and P (even if they do not own them)! */
1515  hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1516  hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1517 
1518  if (P_owns_its_col_starts)
1519  {
1520  hypre_ParCSRMatrixSetColStartsOwner(*P, 1);
1521  }
1522  if (Rt_owns_its_col_starts)
1523  {
1524  hypre_ParCSRMatrixSetColStartsOwner(*Rt, 1);
1525  }
1526 
1527  return new HypreParMatrix(rap);
1528 }
1529 
1531  const Array<int> &ess_dof_list,
1532  const Vector &X, Vector &B)
1533 {
1534  // B -= Ae*X
1535  Ae.Mult(-1.0, X, 1.0, B);
1536 
1537  hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag((hypre_ParCSRMatrix *)A);
1538  double *data = hypre_CSRMatrixData(A_diag);
1539  HYPRE_Int *I = hypre_CSRMatrixI(A_diag);
1540 #ifdef MFEM_DEBUG
1541  HYPRE_Int *J = hypre_CSRMatrixJ(A_diag);
1542  hypre_CSRMatrix *A_offd = hypre_ParCSRMatrixOffd((hypre_ParCSRMatrix *)A);
1543  HYPRE_Int *I_offd = hypre_CSRMatrixI(A_offd);
1544  double *data_offd = hypre_CSRMatrixData(A_offd);
1545 #endif
1546 
1547  for (int i = 0; i < ess_dof_list.Size(); i++)
1548  {
1549  int r = ess_dof_list[i];
1550  B(r) = data[I[r]] * X(r);
1551 #ifdef MFEM_DEBUG
1552  // Check that in the rows specified by the ess_dof_list, the matrix A has
1553  // only one entry -- the diagonal.
1554  // if (I[r+1] != I[r]+1 || J[I[r]] != r || I_offd[r] != I_offd[r+1])
1555  if (J[I[r]] != r)
1556  {
1557  MFEM_ABORT("the diagonal entry must be the first entry in the row!");
1558  }
1559  for (int j = I[r]+1; j < I[r+1]; j++)
1560  {
1561  if (data[j] != 0.0)
1562  {
1563  MFEM_ABORT("all off-diagonal entries must be zero!");
1564  }
1565  }
1566  for (int j = I_offd[r]; j < I_offd[r+1]; j++)
1567  {
1568  if (data_offd[j] != 0.0)
1569  {
1570  MFEM_ABORT("all off-diagonal entries must be zero!");
1571  }
1572  }
1573 #endif
1574  }
1575 }
1576 
1577 // Taubin or "lambda-mu" scheme, which alternates between positive and
1578 // negative step sizes to approximate low-pass filter effect.
1579 
1580 int ParCSRRelax_Taubin(hypre_ParCSRMatrix *A, // matrix to relax with
1581  hypre_ParVector *f, // right-hand side
1582  double lambda,
1583  double mu,
1584  int N,
1585  double max_eig,
1586  hypre_ParVector *u, // initial/updated approximation
1587  hypre_ParVector *r // another temp vector
1588  )
1589 {
1590  hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1591  HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1592 
1593  double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1594  double *r_data = hypre_VectorData(hypre_ParVectorLocalVector(r));
1595 
1596  for (int i = 0; i < N; i++)
1597  {
1598  // get residual: r = f - A*u
1599  hypre_ParVectorCopy(f, r);
1600  hypre_ParCSRMatrixMatvec(-1.0, A, u, 1.0, r);
1601 
1602  double coef;
1603  (0 == (i % 2)) ? coef = lambda : coef = mu;
1604 
1605  for (HYPRE_Int j = 0; j < num_rows; j++)
1606  {
1607  u_data[j] += coef*r_data[j] / max_eig;
1608  }
1609  }
1610 
1611  return 0;
1612 }
1613 
1614 // FIR scheme, which uses Chebyshev polynomials and a window function
1615 // to approximate a low-pass step filter.
1616 
1617 int ParCSRRelax_FIR(hypre_ParCSRMatrix *A, // matrix to relax with
1618  hypre_ParVector *f, // right-hand side
1619  double max_eig,
1620  int poly_order,
1621  double* fir_coeffs,
1622  hypre_ParVector *u, // initial/updated approximation
1623  hypre_ParVector *x0, // temporaries
1624  hypre_ParVector *x1,
1625  hypre_ParVector *x2,
1626  hypre_ParVector *x3)
1627 
1628 {
1629  hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1630  HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1631 
1632  double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1633 
1634  double *x0_data = hypre_VectorData(hypre_ParVectorLocalVector(x0));
1635  double *x1_data = hypre_VectorData(hypre_ParVectorLocalVector(x1));
1636  double *x2_data = hypre_VectorData(hypre_ParVectorLocalVector(x2));
1637  double *x3_data = hypre_VectorData(hypre_ParVectorLocalVector(x3));
1638 
1639  hypre_ParVectorCopy(u, x0);
1640 
1641  // x1 = f -A*x0/max_eig
1642  hypre_ParVectorCopy(f, x1);
1643  hypre_ParCSRMatrixMatvec(-1.0, A, x0, 1.0, x1);
1644 
1645  for (HYPRE_Int i = 0; i < num_rows; i++)
1646  {
1647  x1_data[i] /= -max_eig;
1648  }
1649 
1650  // x1 = x0 -x1
1651  for (HYPRE_Int i = 0; i < num_rows; i++)
1652  {
1653  x1_data[i] = x0_data[i] -x1_data[i];
1654  }
1655 
1656  // x3 = f0*x0 +f1*x1
1657  for (HYPRE_Int i = 0; i < num_rows; i++)
1658  {
1659  x3_data[i] = fir_coeffs[0]*x0_data[i] +fir_coeffs[1]*x1_data[i];
1660  }
1661 
1662  for (int n = 2; n <= poly_order; n++)
1663  {
1664  // x2 = f - A*x1/max_eig
1665  hypre_ParVectorCopy(f, x2);
1666  hypre_ParCSRMatrixMatvec(-1.0, A, x1, 1.0, x2);
1667 
1668  for (HYPRE_Int i = 0; i < num_rows; i++)
1669  {
1670  x2_data[i] /= -max_eig;
1671  }
1672 
1673  // x2 = (x1-x0) +(x1-2*x2)
1674  // x3 = x3 +f[n]*x2
1675  // x0 = x1
1676  // x1 = x2
1677 
1678  for (HYPRE_Int i = 0; i < num_rows; i++)
1679  {
1680  x2_data[i] = (x1_data[i]-x0_data[i]) +(x1_data[i]-2*x2_data[i]);
1681  x3_data[i] += fir_coeffs[n]*x2_data[i];
1682  x0_data[i] = x1_data[i];
1683  x1_data[i] = x2_data[i];
1684  }
1685  }
1686 
1687  for (HYPRE_Int i = 0; i < num_rows; i++)
1688  {
1689  u_data[i] = x3_data[i];
1690  }
1691 
1692  return 0;
1693 }
1694 
1696 {
1697  type = 2;
1698  relax_times = 1;
1699  relax_weight = 1.0;
1700  omega = 1.0;
1701  poly_order = 2;
1702  poly_fraction = .3;
1703  lambda = 0.5;
1704  mu = -0.5;
1705  taubin_iter = 40;
1706 
1707  l1_norms = NULL;
1708  pos_l1_norms = false;
1709  B = X = V = Z = NULL;
1710  X0 = X1 = NULL;
1711  fir_coeffs = NULL;
1712 }
1713 
1715  int _relax_times, double _relax_weight, double _omega,
1716  int _poly_order, double _poly_fraction)
1717 {
1718  type = _type;
1719  relax_times = _relax_times;
1720  relax_weight = _relax_weight;
1721  omega = _omega;
1722  poly_order = _poly_order;
1723  poly_fraction = _poly_fraction;
1724 
1725  l1_norms = NULL;
1726  pos_l1_norms = false;
1727  B = X = V = Z = NULL;
1728  X0 = X1 = NULL;
1729  fir_coeffs = NULL;
1730 
1731  SetOperator(_A);
1732 }
1733 
1734 void HypreSmoother::SetType(HypreSmoother::Type _type, int _relax_times)
1735 {
1736  type = static_cast<int>(_type);
1737  relax_times = _relax_times;
1738 }
1739 
1740 void HypreSmoother::SetSOROptions(double _relax_weight, double _omega)
1741 {
1742  relax_weight = _relax_weight;
1743  omega = _omega;
1744 }
1745 
1746 void HypreSmoother::SetPolyOptions(int _poly_order, double _poly_fraction)
1747 {
1748  poly_order = _poly_order;
1749  poly_fraction = _poly_fraction;
1750 }
1751 
1752 void HypreSmoother::SetTaubinOptions(double _lambda, double _mu,
1753  int _taubin_iter)
1754 {
1755  lambda = _lambda;
1756  mu = _mu;
1757  taubin_iter = _taubin_iter;
1758 }
1759 
1760 void HypreSmoother::SetWindowByName(const char* name)
1761 {
1762  double a = -1, b, c;
1763  if (!strcmp(name,"Rectangular")) { a = 1.0, b = 0.0, c = 0.0; }
1764  if (!strcmp(name,"Hanning")) { a = 0.5, b = 0.5, c = 0.0; }
1765  if (!strcmp(name,"Hamming")) { a = 0.54, b = 0.46, c = 0.0; }
1766  if (!strcmp(name,"Blackman")) { a = 0.42, b = 0.50, c = 0.08; }
1767  if (a < 0)
1768  {
1769  mfem_error("HypreSmoother::SetWindowByName : name not recognized!");
1770  }
1771 
1772  SetWindowParameters(a, b, c);
1773 }
1774 
1775 void HypreSmoother::SetWindowParameters(double a, double b, double c)
1776 {
1777  window_params[0] = a;
1778  window_params[1] = b;
1779  window_params[2] = c;
1780 }
1781 
1783 {
1784  A = const_cast<HypreParMatrix *>(dynamic_cast<const HypreParMatrix *>(&op));
1785  if (A == NULL)
1786  {
1787  mfem_error("HypreSmoother::SetOperator : not HypreParMatrix!");
1788  }
1789 
1790  height = A->Height();
1791  width = A->Width();
1792 
1793  if (B) { delete B; }
1794  if (X) { delete X; }
1795  if (V) { delete V; }
1796  if (Z) { delete Z; }
1797  if (l1_norms)
1798  {
1799  hypre_TFree(l1_norms);
1800  }
1801  delete X0;
1802  delete X1;
1803 
1804  X1 = X0 = Z = V = B = X = NULL;
1805 
1806  if (type >= 1 && type <= 4)
1807  {
1808  hypre_ParCSRComputeL1Norms(*A, type, NULL, &l1_norms);
1809  }
1810  else if (type == 5)
1811  {
1812  l1_norms = hypre_CTAlloc(double, height);
1813  Vector ones(height), diag(l1_norms, height);
1814  ones = 1.0;
1815  A->Mult(ones, diag);
1816  type = 1;
1817  }
1818  else
1819  {
1820  l1_norms = NULL;
1821  }
1822  if (l1_norms && pos_l1_norms)
1823  {
1824  for (int i = 0; i < height; i++)
1825  {
1826  l1_norms[i] = std::abs(l1_norms[i]);
1827  }
1828  }
1829 
1830  if (type == 16)
1831  {
1832  poly_scale = 1;
1833  hypre_ParCSRMaxEigEstimateCG(*A, poly_scale, 10,
1835  Z = new HypreParVector(*A);
1836  }
1837  else if (type == 1001 || type == 1002)
1838  {
1839  poly_scale = 0;
1840  hypre_ParCSRMaxEigEstimateCG(*A, poly_scale, 10,
1842 
1843  // The Taubin and FIR polynomials are defined on [0, 2]
1844  max_eig_est /= 2;
1845 
1846  // Compute window function, Chebyshev coefficients, and allocate temps.
1847  if (type == 1002)
1848  {
1849  // Temporaries for Chebyshev recursive evaluation
1850  Z = new HypreParVector(*A);
1851  X0 = new HypreParVector(*A);
1852  X1 = new HypreParVector(*A);
1853 
1855  }
1856  }
1857 }
1858 
1860 {
1861  if (fir_coeffs)
1862  {
1863  delete [] fir_coeffs;
1864  }
1865 
1866  fir_coeffs = new double[poly_order+1];
1867 
1868  double* window_coeffs = new double[poly_order+1];
1869  double* cheby_coeffs = new double[poly_order+1];
1870 
1871  double a = window_params[0];
1872  double b = window_params[1];
1873  double c = window_params[2];
1874  for (int i = 0; i <= poly_order; i++)
1875  {
1876  double t = (i*M_PI)/(poly_order+1);
1877  window_coeffs[i] = a + b*cos(t) +c*cos(2*t);
1878  }
1879 
1880  double k_pb = poly_fraction*max_eig;
1881  double theta_pb = acos(1.0 -0.5*k_pb);
1882  double sigma = 0.0;
1883  cheby_coeffs[0] = (theta_pb +sigma)/M_PI;
1884  for (int i = 1; i <= poly_order; i++)
1885  {
1886  double t = i*(theta_pb+sigma);
1887  cheby_coeffs[i] = 2.0*sin(t)/(i*M_PI);
1888  }
1889 
1890  for (int i = 0; i <= poly_order; i++)
1891  {
1892  fir_coeffs[i] = window_coeffs[i]*cheby_coeffs[i];
1893  }
1894 
1895  delete[] window_coeffs;
1896  delete[] cheby_coeffs;
1897 }
1898 
1900 {
1901  if (A == NULL)
1902  {
1903  mfem_error("HypreSmoother::Mult (...) : HypreParMatrix A is missing");
1904  return;
1905  }
1906 
1907  if (!iterative_mode)
1908  {
1909  if (type == 0 && relax_times == 1)
1910  {
1911  HYPRE_ParCSRDiagScale(NULL, *A, b, x);
1912  if (relax_weight != 1.0)
1913  {
1914  x *= relax_weight;
1915  }
1916  return;
1917  }
1918  x = 0.0;
1919  }
1920 
1921  if (V == NULL)
1922  {
1923  V = new HypreParVector(*A);
1924  }
1925 
1926  if (type == 1001)
1927  {
1928  for (int sweep = 0; sweep < relax_times; sweep++)
1929  {
1931  max_eig_est,
1932  x, *V);
1933  }
1934  }
1935  else if (type == 1002)
1936  {
1937  for (int sweep = 0; sweep < relax_times; sweep++)
1938  {
1939  ParCSRRelax_FIR(*A, b,
1940  max_eig_est,
1941  poly_order,
1942  fir_coeffs,
1943  x,
1944  *X0, *X1, *V, *Z);
1945  }
1946  }
1947  else
1948  {
1949  if (Z == NULL)
1950  hypre_ParCSRRelax(*A, b, type,
1953  x, *V, NULL);
1954  else
1955  hypre_ParCSRRelax(*A, b, type,
1958  x, *V, *Z);
1959  }
1960 }
1961 
1962 void HypreSmoother::Mult(const Vector &b, Vector &x) const
1963 {
1964  if (A == NULL)
1965  {
1966  mfem_error("HypreSmoother::Mult (...) : HypreParMatrix A is missing");
1967  return;
1968  }
1969  if (B == NULL)
1970  {
1971  B = new HypreParVector(A->GetComm(),
1972  A -> GetGlobalNumRows(),
1973  b.GetData(),
1974  A -> GetRowStarts());
1975  X = new HypreParVector(A->GetComm(),
1976  A -> GetGlobalNumCols(),
1977  x.GetData(),
1978  A -> GetColStarts());
1979  }
1980  else
1981  {
1982  B -> SetData(b.GetData());
1983  X -> SetData(x.GetData());
1984  }
1985 
1986  Mult(*B, *X);
1987 }
1988 
1990 {
1991  if (B) { delete B; }
1992  if (X) { delete X; }
1993  if (V) { delete V; }
1994  if (Z) { delete Z; }
1995  if (l1_norms)
1996  {
1997  hypre_TFree(l1_norms);
1998  }
1999  if (fir_coeffs)
2000  {
2001  delete [] fir_coeffs;
2002  }
2003  if (X0) { delete X0; }
2004  if (X1) { delete X1; }
2005 }
2006 
2007 
2009 {
2010  A = NULL;
2011  setup_called = 0;
2012  B = X = NULL;
2013 }
2014 
2016  : Solver(_A->Height(), _A->Width())
2017 {
2018  A = _A;
2019  setup_called = 0;
2020  B = X = NULL;
2021 }
2022 
2024 {
2025  if (A == NULL)
2026  {
2027  mfem_error("HypreSolver::Mult (...) : HypreParMatrix A is missing");
2028  return;
2029  }
2030  if (!setup_called)
2031  {
2032  SetupFcn()(*this, *A, b, x);
2033  setup_called = 1;
2034  }
2035 
2036  if (!iterative_mode)
2037  {
2038  x = 0.0;
2039  }
2040  SolveFcn()(*this, *A, b, x);
2041 }
2042 
2043 void HypreSolver::Mult(const Vector &b, Vector &x) const
2044 {
2045  if (A == NULL)
2046  {
2047  mfem_error("HypreSolver::Mult (...) : HypreParMatrix A is missing");
2048  return;
2049  }
2050  if (B == NULL)
2051  {
2052  B = new HypreParVector(A->GetComm(),
2053  A -> GetGlobalNumRows(),
2054  b.GetData(),
2055  A -> GetRowStarts());
2056  X = new HypreParVector(A->GetComm(),
2057  A -> GetGlobalNumCols(),
2058  x.GetData(),
2059  A -> GetColStarts());
2060  }
2061  else
2062  {
2063  B -> SetData(b.GetData());
2064  X -> SetData(x.GetData());
2065  }
2066 
2067  Mult(*B, *X);
2068 }
2069 
2071 {
2072  if (B) { delete B; }
2073  if (X) { delete X; }
2074 }
2075 
2076 
2078 {
2079  MPI_Comm comm;
2080 
2081  iterative_mode = true;
2082 
2083  HYPRE_ParCSRMatrixGetComm(*A, &comm);
2084 
2085  HYPRE_ParCSRPCGCreate(comm, &pcg_solver);
2086 }
2087 
2088 void HyprePCG::SetTol(double tol)
2089 {
2090  HYPRE_PCGSetTol(pcg_solver, tol);
2091 }
2092 
2093 void HyprePCG::SetMaxIter(int max_iter)
2094 {
2095  HYPRE_PCGSetMaxIter(pcg_solver, max_iter);
2096 }
2097 
2098 void HyprePCG::SetLogging(int logging)
2099 {
2100  HYPRE_PCGSetLogging(pcg_solver, logging);
2101 }
2102 
2103 void HyprePCG::SetPrintLevel(int print_lvl)
2104 {
2105  HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_lvl);
2106 }
2107 
2109 {
2110  HYPRE_ParCSRPCGSetPrecond(pcg_solver,
2111  precond.SolveFcn(),
2112  precond.SetupFcn(),
2113  precond);
2114 }
2115 
2116 void HyprePCG::SetResidualConvergenceOptions(int res_frequency, double rtol)
2117 {
2118  HYPRE_PCGSetTwoNorm(pcg_solver, 1);
2119  if (res_frequency > 0)
2120  {
2121  HYPRE_PCGSetRecomputeResidualP(pcg_solver, res_frequency);
2122  }
2123  if (rtol > 0.0)
2124  {
2125  HYPRE_PCGSetResidualTol(pcg_solver, rtol);
2126  }
2127 }
2128 
2130 {
2131  int myid;
2132  HYPRE_Int time_index = 0;
2133  HYPRE_Int num_iterations;
2134  double final_res_norm;
2135  MPI_Comm comm;
2136  HYPRE_Int print_level;
2137 
2138  HYPRE_PCGGetPrintLevel(pcg_solver, &print_level);
2139  HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level%3);
2140 
2141  HYPRE_ParCSRMatrixGetComm(*A, &comm);
2142 
2143  if (!setup_called)
2144  {
2145  if (print_level > 0 && print_level < 3)
2146  {
2147  time_index = hypre_InitializeTiming("PCG Setup");
2148  hypre_BeginTiming(time_index);
2149  }
2150 
2151  HYPRE_ParCSRPCGSetup(pcg_solver, *A, b, x);
2152  setup_called = 1;
2153 
2154  if (print_level > 0 && print_level < 3)
2155  {
2156  hypre_EndTiming(time_index);
2157  hypre_PrintTiming("Setup phase times", comm);
2158  hypre_FinalizeTiming(time_index);
2159  hypre_ClearTiming();
2160  }
2161  }
2162 
2163  if (print_level > 0 && print_level < 3)
2164  {
2165  time_index = hypre_InitializeTiming("PCG Solve");
2166  hypre_BeginTiming(time_index);
2167  }
2168 
2169  if (!iterative_mode)
2170  {
2171  x = 0.0;
2172  }
2173 
2174  HYPRE_ParCSRPCGSolve(pcg_solver, *A, b, x);
2175 
2176  if (print_level > 0)
2177  {
2178  if (print_level < 3)
2179  {
2180  hypre_EndTiming(time_index);
2181  hypre_PrintTiming("Solve phase times", comm);
2182  hypre_FinalizeTiming(time_index);
2183  hypre_ClearTiming();
2184  }
2185 
2186  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_iterations);
2187  HYPRE_ParCSRPCGGetFinalRelativeResidualNorm(pcg_solver,
2188  &final_res_norm);
2189 
2190  MPI_Comm_rank(comm, &myid);
2191 
2192  if (myid == 0)
2193  {
2194  mfem::out << "PCG Iterations = " << num_iterations << endl
2195  << "Final PCG Relative Residual Norm = " << final_res_norm
2196  << endl;
2197  }
2198  }
2199  HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level);
2200 }
2201 
2203 {
2204  HYPRE_ParCSRPCGDestroy(pcg_solver);
2205 }
2206 
2207 
2209 {
2210  MPI_Comm comm;
2211 
2212  int k_dim = 50;
2213  int max_iter = 100;
2214  double tol = 1e-6;
2215 
2216  iterative_mode = true;
2217 
2218  HYPRE_ParCSRMatrixGetComm(*A, &comm);
2219 
2220  HYPRE_ParCSRGMRESCreate(comm, &gmres_solver);
2221  HYPRE_ParCSRGMRESSetKDim(gmres_solver, k_dim);
2222  HYPRE_ParCSRGMRESSetMaxIter(gmres_solver, max_iter);
2223  HYPRE_ParCSRGMRESSetTol(gmres_solver, tol);
2224 }
2225 
2226 void HypreGMRES::SetTol(double tol)
2227 {
2228  HYPRE_GMRESSetTol(gmres_solver, tol);
2229 }
2230 
2231 void HypreGMRES::SetMaxIter(int max_iter)
2232 {
2233  HYPRE_GMRESSetMaxIter(gmres_solver, max_iter);
2234 }
2235 
2236 void HypreGMRES::SetKDim(int k_dim)
2237 {
2238  HYPRE_GMRESSetKDim(gmres_solver, k_dim);
2239 }
2240 
2241 void HypreGMRES::SetLogging(int logging)
2242 {
2243  HYPRE_GMRESSetLogging(gmres_solver, logging);
2244 }
2245 
2246 void HypreGMRES::SetPrintLevel(int print_lvl)
2247 {
2248  HYPRE_GMRESSetPrintLevel(gmres_solver, print_lvl);
2249 }
2250 
2252 {
2253  HYPRE_ParCSRGMRESSetPrecond(gmres_solver,
2254  precond.SolveFcn(),
2255  precond.SetupFcn(),
2256  precond);
2257 }
2258 
2260 {
2261  int myid;
2262  HYPRE_Int time_index = 0;
2263  HYPRE_Int num_iterations;
2264  double final_res_norm;
2265  MPI_Comm comm;
2266  HYPRE_Int print_level;
2267 
2268  HYPRE_GMRESGetPrintLevel(gmres_solver, &print_level);
2269 
2270  HYPRE_ParCSRMatrixGetComm(*A, &comm);
2271 
2272  if (!setup_called)
2273  {
2274  if (print_level > 0)
2275  {
2276  time_index = hypre_InitializeTiming("GMRES Setup");
2277  hypre_BeginTiming(time_index);
2278  }
2279 
2280  HYPRE_ParCSRGMRESSetup(gmres_solver, *A, b, x);
2281  setup_called = 1;
2282 
2283  if (print_level > 0)
2284  {
2285  hypre_EndTiming(time_index);
2286  hypre_PrintTiming("Setup phase times", comm);
2287  hypre_FinalizeTiming(time_index);
2288  hypre_ClearTiming();
2289  }
2290  }
2291 
2292  if (print_level > 0)
2293  {
2294  time_index = hypre_InitializeTiming("GMRES Solve");
2295  hypre_BeginTiming(time_index);
2296  }
2297 
2298  if (!iterative_mode)
2299  {
2300  x = 0.0;
2301  }
2302 
2303  HYPRE_ParCSRGMRESSolve(gmres_solver, *A, b, x);
2304 
2305  if (print_level > 0)
2306  {
2307  hypre_EndTiming(time_index);
2308  hypre_PrintTiming("Solve phase times", comm);
2309  hypre_FinalizeTiming(time_index);
2310  hypre_ClearTiming();
2311 
2312  HYPRE_ParCSRGMRESGetNumIterations(gmres_solver, &num_iterations);
2313  HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(gmres_solver,
2314  &final_res_norm);
2315 
2316  MPI_Comm_rank(comm, &myid);
2317 
2318  if (myid == 0)
2319  {
2320  mfem::out << "GMRES Iterations = " << num_iterations << endl
2321  << "Final GMRES Relative Residual Norm = " << final_res_norm
2322  << endl;
2323  }
2324  }
2325 }
2326 
2328 {
2329  HYPRE_ParCSRGMRESDestroy(gmres_solver);
2330 }
2331 
2332 
2334 {
2335  MPI_Comm comm;
2336 
2337  int sai_max_levels = 1;
2338  double sai_threshold = 0.1;
2339  double sai_filter = 0.1;
2340  int sai_sym = 0;
2341  double sai_loadbal = 0.0;
2342  int sai_reuse = 0;
2343  int sai_logging = 1;
2344 
2345  HYPRE_ParCSRMatrixGetComm(A, &comm);
2346 
2347  HYPRE_ParaSailsCreate(comm, &sai_precond);
2348  HYPRE_ParaSailsSetParams(sai_precond, sai_threshold, sai_max_levels);
2349  HYPRE_ParaSailsSetFilter(sai_precond, sai_filter);
2350  HYPRE_ParaSailsSetSym(sai_precond, sai_sym);
2351  HYPRE_ParaSailsSetLoadbal(sai_precond, sai_loadbal);
2352  HYPRE_ParaSailsSetReuse(sai_precond, sai_reuse);
2353  HYPRE_ParaSailsSetLogging(sai_precond, sai_logging);
2354 }
2355 
2357 {
2358  HYPRE_ParaSailsSetSym(sai_precond, sym);
2359 }
2360 
2362 {
2363  HYPRE_ParaSailsDestroy(sai_precond);
2364 }
2365 
2366 
2368 {
2369  HYPRE_BoomerAMGCreate(&amg_precond);
2370  SetDefaultOptions();
2371 }
2372 
2374 {
2375  HYPRE_BoomerAMGCreate(&amg_precond);
2376  SetDefaultOptions();
2377 }
2378 
2379 void HypreBoomerAMG::SetDefaultOptions()
2380 {
2381  // AMG coarsening options:
2382  int coarsen_type = 10; // 10 = HMIS, 8 = PMIS, 6 = Falgout, 0 = CLJP
2383  int agg_levels = 1; // number of aggressive coarsening levels
2384  double theta = 0.25; // strength threshold: 0.25, 0.5, 0.8
2385 
2386  // AMG interpolation options:
2387  int interp_type = 6; // 6 = extended+i, 0 = classical
2388  int Pmax = 4; // max number of elements per row in P
2389 
2390  // AMG relaxation options:
2391  int relax_type = 8; // 8 = l1-GS, 6 = symm. GS, 3 = GS, 18 = l1-Jacobi
2392  int relax_sweeps = 1; // relaxation sweeps on each level
2393 
2394  // Additional options:
2395  int print_level = 1; // print AMG iterations? 1 = no, 2 = yes
2396  int max_levels = 25; // max number of levels in AMG hierarchy
2397 
2398  HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2399  HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2400  HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2401  HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2402  HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2403  HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2404  HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2405  HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2406  HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels);
2407 
2408  // Use as a preconditioner (one V-cycle, zero tolerance)
2409  HYPRE_BoomerAMGSetMaxIter(amg_precond, 1);
2410  HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2411 }
2412 
2413 void HypreBoomerAMG::ResetAMGPrecond()
2414 {
2415  HYPRE_Int coarsen_type;
2416  HYPRE_Int agg_levels;
2417  HYPRE_Int relax_type;
2418  HYPRE_Int relax_sweeps;
2419  HYPRE_Real theta;
2420  HYPRE_Int interp_type;
2421  HYPRE_Int Pmax;
2422  HYPRE_Int print_level;
2423  HYPRE_Int dim;
2424  HYPRE_Int nrbms = rbms.Size();
2425  HYPRE_Int nodal;
2426  HYPRE_Int nodal_diag;
2427  HYPRE_Int relax_coarse;
2428  HYPRE_Int interp_vec_variant;
2429  HYPRE_Int q_max;
2430  HYPRE_Int smooth_interp_vectors;
2431  HYPRE_Int interp_refine;
2432 
2433  hypre_ParAMGData *amg_data = (hypre_ParAMGData *)amg_precond;
2434 
2435  // read options from amg_precond
2436  HYPRE_BoomerAMGGetCoarsenType(amg_precond, &coarsen_type);
2437  agg_levels = hypre_ParAMGDataAggNumLevels(amg_data);
2438  relax_type = hypre_ParAMGDataUserRelaxType(amg_data);
2439  relax_sweeps = hypre_ParAMGDataUserNumSweeps(amg_data);
2440  HYPRE_BoomerAMGGetStrongThreshold(amg_precond, &theta);
2441  hypre_BoomerAMGGetInterpType(amg_precond, &interp_type);
2442  HYPRE_BoomerAMGGetPMaxElmts(amg_precond, &Pmax);
2443  HYPRE_BoomerAMGGetPrintLevel(amg_precond, &print_level);
2444  HYPRE_BoomerAMGGetNumFunctions(amg_precond, &dim);
2445  if (nrbms) // elasticity solver options
2446  {
2447  nodal = hypre_ParAMGDataNodal(amg_data);
2448  nodal_diag = hypre_ParAMGDataNodalDiag(amg_data);
2449  HYPRE_BoomerAMGGetCycleRelaxType(amg_precond, &relax_coarse, 3);
2450  interp_vec_variant = hypre_ParAMGInterpVecVariant(amg_data);
2451  q_max = hypre_ParAMGInterpVecQMax(amg_data);
2452  smooth_interp_vectors = hypre_ParAMGSmoothInterpVectors(amg_data);
2453  interp_refine = hypre_ParAMGInterpRefine(amg_data);
2454  }
2455 
2456  HYPRE_BoomerAMGDestroy(amg_precond);
2457  HYPRE_BoomerAMGCreate(&amg_precond);
2458 
2459  HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2460  HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2461  HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2462  HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2463  HYPRE_BoomerAMGSetMaxLevels(amg_precond, 25);
2464  HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2465  HYPRE_BoomerAMGSetMaxIter(amg_precond, 1); // one V-cycle
2466  HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2467  HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2468  HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2469  HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2470  HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2471  if (nrbms)
2472  {
2473  HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2474  HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2475  HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2476  HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2477  HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2478  HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2479  HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2480  RecomputeRBMs();
2481  HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.Size(), rbms.GetData());
2482  }
2483 }
2484 
2486 {
2487  const HypreParMatrix *new_A = dynamic_cast<const HypreParMatrix *>(&op);
2488  MFEM_VERIFY(new_A, "new Operator must be a HypreParMatrix!");
2489 
2490  if (A) { ResetAMGPrecond(); }
2491 
2492  // update base classes: Operator, Solver, HypreSolver
2493  height = new_A->Height();
2494  width = new_A->Width();
2495  A = const_cast<HypreParMatrix *>(new_A);
2496  setup_called = 0;
2497  delete X;
2498  delete B;
2499  B = X = NULL;
2500 }
2501 
2503 {
2504  HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2505 
2506  // More robust options with respect to convergence
2507  HYPRE_BoomerAMGSetAggNumLevels(amg_precond, 0);
2508  HYPRE_BoomerAMGSetStrongThreshold(amg_precond, 0.5);
2509 }
2510 
2511 // Rotational rigid-body mode functions, used in SetElasticityOptions()
2512 static void func_rxy(const Vector &x, Vector &y)
2513 {
2514  y = 0.0; y(0) = x(1); y(1) = -x(0);
2515 }
2516 static void func_ryz(const Vector &x, Vector &y)
2517 {
2518  y = 0.0; y(1) = x(2); y(2) = -x(1);
2519 }
2520 static void func_rzx(const Vector &x, Vector &y)
2521 {
2522  y = 0.0; y(2) = x(0); y(0) = -x(2);
2523 }
2524 
2525 void HypreBoomerAMG::RecomputeRBMs()
2526 {
2527  int nrbms;
2528  Array<HypreParVector*> gf_rbms;
2529  int dim = fespace->GetParMesh()->Dimension();
2530 
2531  for (int i = 0; i < rbms.Size(); i++)
2532  {
2533  HYPRE_ParVectorDestroy(rbms[i]);
2534  }
2535 
2536  if (dim == 2)
2537  {
2538  nrbms = 1;
2539 
2540  VectorFunctionCoefficient coeff_rxy(2, func_rxy);
2541 
2542  ParGridFunction rbms_rxy(fespace);
2543  rbms_rxy.ProjectCoefficient(coeff_rxy);
2544 
2545  rbms.SetSize(nrbms);
2546  gf_rbms.SetSize(nrbms);
2547  gf_rbms[0] = rbms_rxy.ParallelAverage();
2548  }
2549  else if (dim == 3)
2550  {
2551  nrbms = 3;
2552 
2553  VectorFunctionCoefficient coeff_rxy(3, func_rxy);
2554  VectorFunctionCoefficient coeff_ryz(3, func_ryz);
2555  VectorFunctionCoefficient coeff_rzx(3, func_rzx);
2556 
2557  ParGridFunction rbms_rxy(fespace);
2558  ParGridFunction rbms_ryz(fespace);
2559  ParGridFunction rbms_rzx(fespace);
2560  rbms_rxy.ProjectCoefficient(coeff_rxy);
2561  rbms_ryz.ProjectCoefficient(coeff_ryz);
2562  rbms_rzx.ProjectCoefficient(coeff_rzx);
2563 
2564  rbms.SetSize(nrbms);
2565  gf_rbms.SetSize(nrbms);
2566  gf_rbms[0] = rbms_rxy.ParallelAverage();
2567  gf_rbms[1] = rbms_ryz.ParallelAverage();
2568  gf_rbms[2] = rbms_rzx.ParallelAverage();
2569  }
2570  else
2571  {
2572  nrbms = 0;
2573  rbms.SetSize(nrbms);
2574  }
2575 
2576  // Transfer the RBMs from the ParGridFunction to the HYPRE_ParVector objects
2577  for (int i = 0; i < nrbms; i++)
2578  {
2579  rbms[i] = gf_rbms[i]->StealParVector();
2580  delete gf_rbms[i];
2581  }
2582 }
2583 
2585 {
2586  // Save the finite element space to support multiple calls to SetOperator()
2587  this->fespace = fespace;
2588 
2589  // Make sure the systems AMG options are set
2590  int dim = fespace->GetParMesh()->Dimension();
2591  SetSystemsOptions(dim);
2592 
2593  // Nodal coarsening options (nodal coarsening is required for this solver)
2594  // See hypre's new_ij driver and the paper for descriptions.
2595  int nodal = 4; // strength reduction norm: 1, 3 or 4
2596  int nodal_diag = 1; // diagonal in strength matrix: 0, 1 or 2
2597  int relax_coarse = 8; // smoother on the coarsest grid: 8, 99 or 29
2598 
2599  // Elasticity interpolation options
2600  int interp_vec_variant = 2; // 1 = GM-1, 2 = GM-2, 3 = LN
2601  int q_max = 4; // max elements per row for each Q
2602  int smooth_interp_vectors = 1; // smooth the rigid-body modes?
2603 
2604  // Optionally pre-process the interpolation matrix through iterative weight
2605  // refinement (this is generally applicable for any system)
2606  int interp_refine = 1;
2607 
2608  HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2609  HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2610  HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2611  HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2612  HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2613  HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2614  HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2615 
2616  RecomputeRBMs();
2617  HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.Size(), rbms.GetData());
2618 }
2619 
2621 {
2622  for (int i = 0; i < rbms.Size(); i++)
2623  {
2624  HYPRE_ParVectorDestroy(rbms[i]);
2625  }
2626 
2627  HYPRE_BoomerAMGDestroy(amg_precond);
2628 }
2629 
2630 
2632  : HypreSolver(&A)
2633 {
2634  int cycle_type = 13;
2635  int rlx_type = 2;
2636  int rlx_sweeps = 1;
2637  double rlx_weight = 1.0;
2638  double rlx_omega = 1.0;
2639  int amg_coarsen_type = 10;
2640  int amg_agg_levels = 1;
2641  int amg_rlx_type = 8;
2642  double theta = 0.25;
2643  int amg_interp_type = 6;
2644  int amg_Pmax = 4;
2645 
2646  int dim = edge_fespace->GetMesh()->Dimension();
2647  int sdim = edge_fespace->GetMesh()->SpaceDimension();
2648  const FiniteElementCollection *edge_fec = edge_fespace->FEColl();
2649 
2650  bool trace_space, rt_trace_space;
2651  ND_Trace_FECollection *nd_tr_fec = NULL;
2652  trace_space = dynamic_cast<const ND_Trace_FECollection*>(edge_fec);
2653  rt_trace_space = dynamic_cast<const RT_Trace_FECollection*>(edge_fec);
2654  trace_space = trace_space || rt_trace_space;
2655 
2656  int p = 1;
2657  if (edge_fespace->GetNE() > 0)
2658  {
2659  if (trace_space)
2660  {
2661  p = edge_fespace->GetFaceOrder(0);
2662  if (dim == 2) { p++; }
2663  }
2664  else
2665  {
2666  p = edge_fespace->GetOrder(0);
2667  }
2668  }
2669 
2670  ParMesh *pmesh = edge_fespace->GetParMesh();
2671  if (rt_trace_space)
2672  {
2673  nd_tr_fec = new ND_Trace_FECollection(p, dim);
2674  edge_fespace = new ParFiniteElementSpace(pmesh, nd_tr_fec);
2675  }
2676 
2677  HYPRE_AMSCreate(&ams);
2678 
2679  HYPRE_AMSSetDimension(ams, sdim); // 2D H(div) and 3D H(curl) problems
2680  HYPRE_AMSSetTol(ams, 0.0);
2681  HYPRE_AMSSetMaxIter(ams, 1); // use as a preconditioner
2682  HYPRE_AMSSetCycleType(ams, cycle_type);
2683  HYPRE_AMSSetPrintLevel(ams, 1);
2684 
2685  // define the nodal linear finite element space associated with edge_fespace
2686  FiniteElementCollection *vert_fec;
2687  if (trace_space)
2688  {
2689  vert_fec = new H1_Trace_FECollection(p, dim);
2690  }
2691  else
2692  {
2693  vert_fec = new H1_FECollection(p, dim);
2694  }
2695  ParFiniteElementSpace *vert_fespace = new ParFiniteElementSpace(pmesh,
2696  vert_fec);
2697 
2698  // generate and set the vertex coordinates
2699  if (p == 1)
2700  {
2701  ParGridFunction x_coord(vert_fespace);
2702  ParGridFunction y_coord(vert_fespace);
2703  ParGridFunction z_coord(vert_fespace);
2704  double *coord;
2705  for (int i = 0; i < pmesh->GetNV(); i++)
2706  {
2707  coord = pmesh -> GetVertex(i);
2708  x_coord(i) = coord[0];
2709  y_coord(i) = coord[1];
2710  if (sdim == 3) { z_coord(i) = coord[2]; }
2711  }
2712  x = x_coord.ParallelProject();
2713  y = y_coord.ParallelProject();
2714  if (sdim == 2)
2715  {
2716  z = NULL;
2717  HYPRE_AMSSetCoordinateVectors(ams, *x, *y, NULL);
2718  }
2719  else
2720  {
2721  z = z_coord.ParallelProject();
2722  HYPRE_AMSSetCoordinateVectors(ams, *x, *y, *z);
2723  }
2724  }
2725  else
2726  {
2727  x = NULL;
2728  y = NULL;
2729  z = NULL;
2730  }
2731 
2732  // generate and set the discrete gradient
2734  grad = new ParDiscreteLinearOperator(vert_fespace, edge_fespace);
2735  if (trace_space)
2736  {
2738  }
2739  else
2740  {
2742  }
2743  grad->Assemble();
2744  grad->Finalize();
2745  G = grad->ParallelAssemble();
2746  HYPRE_AMSSetDiscreteGradient(ams, *G);
2747  delete grad;
2748 
2749  // generate and set the Nedelec interpolation matrices
2750  Pi = Pix = Piy = Piz = NULL;
2751  if (p > 1)
2752  {
2753  ParFiniteElementSpace *vert_fespace_d
2754  = new ParFiniteElementSpace(pmesh, vert_fec, sdim, Ordering::byVDIM);
2755 
2757  id_ND = new ParDiscreteLinearOperator(vert_fespace_d, edge_fespace);
2758  if (trace_space)
2759  {
2761  }
2762  else
2763  {
2765  }
2766  id_ND->Assemble();
2767  id_ND->Finalize();
2768 
2769  if (cycle_type < 10)
2770  {
2771  Pi = id_ND->ParallelAssemble();
2772  }
2773  else
2774  {
2775  Array2D<HypreParMatrix *> Pi_blocks;
2776  id_ND->GetParBlocks(Pi_blocks);
2777  Pix = Pi_blocks(0,0);
2778  Piy = Pi_blocks(0,1);
2779  if (sdim == 3) { Piz = Pi_blocks(0,2); }
2780  }
2781 
2782  delete id_ND;
2783 
2784  HYPRE_ParCSRMatrix HY_Pi = (Pi) ? (HYPRE_ParCSRMatrix) *Pi : NULL;
2785  HYPRE_ParCSRMatrix HY_Pix = (Pix) ? (HYPRE_ParCSRMatrix) *Pix : NULL;
2786  HYPRE_ParCSRMatrix HY_Piy = (Piy) ? (HYPRE_ParCSRMatrix) *Piy : NULL;
2787  HYPRE_ParCSRMatrix HY_Piz = (Piz) ? (HYPRE_ParCSRMatrix) *Piz : NULL;
2788  HYPRE_AMSSetInterpolations(ams, HY_Pi, HY_Pix, HY_Piy, HY_Piz);
2789 
2790  delete vert_fespace_d;
2791  }
2792 
2793  delete vert_fespace;
2794  delete vert_fec;
2795 
2796  if (rt_trace_space)
2797  {
2798  delete edge_fespace;
2799  delete nd_tr_fec;
2800  }
2801 
2802  // set additional AMS options
2803  HYPRE_AMSSetSmoothingOptions(ams, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
2804  HYPRE_AMSSetAlphaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2805  theta, amg_interp_type, amg_Pmax);
2806  HYPRE_AMSSetBetaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2807  theta, amg_interp_type, amg_Pmax);
2808 }
2809 
2811 {
2812  HYPRE_AMSDestroy(ams);
2813 
2814  delete x;
2815  delete y;
2816  delete z;
2817 
2818  delete G;
2819  delete Pi;
2820  delete Pix;
2821  delete Piy;
2822  delete Piz;
2823 }
2824 
2825 void HypreAMS::SetPrintLevel(int print_lvl)
2826 {
2827  HYPRE_AMSSetPrintLevel(ams, print_lvl);
2828 }
2829 
2831  : HypreSolver(&A)
2832 {
2833  int cycle_type = 11;
2834  int rlx_type = 2;
2835  int rlx_sweeps = 1;
2836  double rlx_weight = 1.0;
2837  double rlx_omega = 1.0;
2838  int amg_coarsen_type = 10;
2839  int amg_agg_levels = 1;
2840  int amg_rlx_type = 8;
2841  double theta = 0.25;
2842  int amg_interp_type = 6;
2843  int amg_Pmax = 4;
2844  int ams_cycle_type = 14;
2845 
2846  const FiniteElementCollection *face_fec = face_fespace->FEColl();
2847  bool trace_space =
2848  (dynamic_cast<const RT_Trace_FECollection*>(face_fec) != NULL);
2849  int p = 1;
2850  if (face_fespace->GetNE() > 0)
2851  {
2852  if (trace_space)
2853  {
2854  p = face_fespace->GetFaceOrder(0) + 1;
2855  }
2856  else
2857  {
2858  p = face_fespace->GetOrder(0);
2859  }
2860  }
2861 
2862  HYPRE_ADSCreate(&ads);
2863 
2864  HYPRE_ADSSetTol(ads, 0.0);
2865  HYPRE_ADSSetMaxIter(ads, 1); // use as a preconditioner
2866  HYPRE_ADSSetCycleType(ads, cycle_type);
2867  HYPRE_ADSSetPrintLevel(ads, 1);
2868 
2869  // define the nodal and edge finite element spaces associated with face_fespace
2870  ParMesh *pmesh = (ParMesh *) face_fespace->GetMesh();
2871  FiniteElementCollection *vert_fec, *edge_fec;
2872  if (trace_space)
2873  {
2874  vert_fec = new H1_Trace_FECollection(p, 3);
2875  edge_fec = new ND_Trace_FECollection(p, 3);
2876  }
2877  else
2878  {
2879  vert_fec = new H1_FECollection(p, 3);
2880  edge_fec = new ND_FECollection(p, 3);
2881  }
2882 
2883  ParFiniteElementSpace *vert_fespace = new ParFiniteElementSpace(pmesh,
2884  vert_fec);
2885  ParFiniteElementSpace *edge_fespace = new ParFiniteElementSpace(pmesh,
2886  edge_fec);
2887 
2888  // generate and set the vertex coordinates
2889  if (p == 1)
2890  {
2891  ParGridFunction x_coord(vert_fespace);
2892  ParGridFunction y_coord(vert_fespace);
2893  ParGridFunction z_coord(vert_fespace);
2894  double *coord;
2895  for (int i = 0; i < pmesh->GetNV(); i++)
2896  {
2897  coord = pmesh -> GetVertex(i);
2898  x_coord(i) = coord[0];
2899  y_coord(i) = coord[1];
2900  z_coord(i) = coord[2];
2901  }
2902  x = x_coord.ParallelProject();
2903  y = y_coord.ParallelProject();
2904  z = z_coord.ParallelProject();
2905  HYPRE_ADSSetCoordinateVectors(ads, *x, *y, *z);
2906  }
2907  else
2908  {
2909  x = NULL;
2910  y = NULL;
2911  z = NULL;
2912  }
2913 
2914  // generate and set the discrete curl
2916  curl = new ParDiscreteLinearOperator(edge_fespace, face_fespace);
2917  if (trace_space)
2918  {
2920  }
2921  else
2922  {
2924  }
2925  curl->Assemble();
2926  curl->Finalize();
2927  C = curl->ParallelAssemble();
2928  C->CopyColStarts(); // since we'll delete edge_fespace
2929  HYPRE_ADSSetDiscreteCurl(ads, *C);
2930  delete curl;
2931 
2932  // generate and set the discrete gradient
2934  grad = new ParDiscreteLinearOperator(vert_fespace, edge_fespace);
2935  if (trace_space)
2936  {
2938  }
2939  else
2940  {
2942  }
2943  grad->Assemble();
2944  grad->Finalize();
2945  G = grad->ParallelAssemble();
2946  G->CopyColStarts(); // since we'll delete vert_fespace
2947  G->CopyRowStarts(); // since we'll delete edge_fespace
2948  HYPRE_ADSSetDiscreteGradient(ads, *G);
2949  delete grad;
2950 
2951  // generate and set the Nedelec and Raviart-Thomas interpolation matrices
2952  RT_Pi = RT_Pix = RT_Piy = RT_Piz = NULL;
2953  ND_Pi = ND_Pix = ND_Piy = ND_Piz = NULL;
2954  if (p > 1)
2955  {
2956  ParFiniteElementSpace *vert_fespace_d
2957  = new ParFiniteElementSpace(pmesh, vert_fec, 3, Ordering::byVDIM);
2958 
2960  id_ND = new ParDiscreteLinearOperator(vert_fespace_d, edge_fespace);
2961  if (trace_space)
2962  {
2964  }
2965  else
2966  {
2968  }
2969  id_ND->Assemble();
2970  id_ND->Finalize();
2971 
2972  if (ams_cycle_type < 10)
2973  {
2974  ND_Pi = id_ND->ParallelAssemble();
2975  ND_Pi->CopyColStarts(); // since we'll delete vert_fespace_d
2976  ND_Pi->CopyRowStarts(); // since we'll delete edge_fespace
2977  }
2978  else
2979  {
2980  Array2D<HypreParMatrix *> ND_Pi_blocks;
2981  id_ND->GetParBlocks(ND_Pi_blocks);
2982  ND_Pix = ND_Pi_blocks(0,0);
2983  ND_Piy = ND_Pi_blocks(0,1);
2984  ND_Piz = ND_Pi_blocks(0,2);
2985  }
2986 
2987  delete id_ND;
2988 
2990  id_RT = new ParDiscreteLinearOperator(vert_fespace_d, face_fespace);
2991  if (trace_space)
2992  {
2994  }
2995  else
2996  {
2998  }
2999  id_RT->Assemble();
3000  id_RT->Finalize();
3001 
3002  if (cycle_type < 10)
3003  {
3004  RT_Pi = id_RT->ParallelAssemble();
3005  RT_Pi->CopyColStarts(); // since we'll delete vert_fespace_d
3006  }
3007  else
3008  {
3009  Array2D<HypreParMatrix *> RT_Pi_blocks;
3010  id_RT->GetParBlocks(RT_Pi_blocks);
3011  RT_Pix = RT_Pi_blocks(0,0);
3012  RT_Piy = RT_Pi_blocks(0,1);
3013  RT_Piz = RT_Pi_blocks(0,2);
3014  }
3015 
3016  delete id_RT;
3017 
3018  HYPRE_ParCSRMatrix HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz;
3019  HY_RT_Pi = (RT_Pi) ? (HYPRE_ParCSRMatrix) *RT_Pi : NULL;
3020  HY_RT_Pix = (RT_Pix) ? (HYPRE_ParCSRMatrix) *RT_Pix : NULL;
3021  HY_RT_Piy = (RT_Piy) ? (HYPRE_ParCSRMatrix) *RT_Piy : NULL;
3022  HY_RT_Piz = (RT_Piz) ? (HYPRE_ParCSRMatrix) *RT_Piz : NULL;
3023  HYPRE_ParCSRMatrix HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz;
3024  HY_ND_Pi = (ND_Pi) ? (HYPRE_ParCSRMatrix) *ND_Pi : NULL;
3025  HY_ND_Pix = (ND_Pix) ? (HYPRE_ParCSRMatrix) *ND_Pix : NULL;
3026  HY_ND_Piy = (ND_Piy) ? (HYPRE_ParCSRMatrix) *ND_Piy : NULL;
3027  HY_ND_Piz = (ND_Piz) ? (HYPRE_ParCSRMatrix) *ND_Piz : NULL;
3028  HYPRE_ADSSetInterpolations(ads,
3029  HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz,
3030  HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz);
3031 
3032  delete vert_fespace_d;
3033  }
3034 
3035  delete vert_fec;
3036  delete vert_fespace;
3037  delete edge_fec;
3038  delete edge_fespace;
3039 
3040  // set additional ADS options
3041  HYPRE_ADSSetSmoothingOptions(ads, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
3042  HYPRE_ADSSetAMGOptions(ads, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
3043  theta, amg_interp_type, amg_Pmax);
3044  HYPRE_ADSSetAMSOptions(ads, ams_cycle_type, amg_coarsen_type, amg_agg_levels,
3045  amg_rlx_type, theta, amg_interp_type, amg_Pmax);
3046 }
3047 
3049 {
3050  HYPRE_ADSDestroy(ads);
3051 
3052  delete x;
3053  delete y;
3054  delete z;
3055 
3056  delete G;
3057  delete C;
3058 
3059  delete RT_Pi;
3060  delete RT_Pix;
3061  delete RT_Piy;
3062  delete RT_Piz;
3063 
3064  delete ND_Pi;
3065  delete ND_Pix;
3066  delete ND_Piy;
3067  delete ND_Piz;
3068 }
3069 
3070 void HypreADS::SetPrintLevel(int print_lvl)
3071 {
3072  HYPRE_ADSSetPrintLevel(ads, print_lvl);
3073 }
3074 
3075 HypreLOBPCG::HypreMultiVector::HypreMultiVector(int n, HypreParVector & v,
3076  mv_InterfaceInterpreter & interpreter)
3077  : hpv(NULL),
3078  nv(n)
3079 {
3080  mv_ptr = mv_MultiVectorCreateFromSampleVector(&interpreter, nv,
3081  (HYPRE_ParVector)v);
3082 
3083  HYPRE_ParVector* vecs = NULL;
3084  {
3085  mv_TempMultiVector* tmp =
3086  (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
3087  vecs = (HYPRE_ParVector*)(tmp -> vector);
3088  }
3089 
3090  hpv = new HypreParVector*[nv];
3091  for (int i=0; i<nv; i++)
3092  {
3093  hpv[i] = new HypreParVector(vecs[i]);
3094  }
3095 }
3096 
3097 HypreLOBPCG::HypreMultiVector::~HypreMultiVector()
3098 {
3099  if ( hpv != NULL )
3100  {
3101  for (int i=0; i<nv; i++)
3102  {
3103  delete hpv[i];
3104  }
3105  delete [] hpv;
3106  }
3107 
3108  mv_MultiVectorDestroy(mv_ptr);
3109 }
3110 
3111 void
3112 HypreLOBPCG::HypreMultiVector::Randomize(HYPRE_Int seed)
3113 {
3114  mv_MultiVectorSetRandom(mv_ptr, seed);
3115 }
3116 
3117 HypreParVector &
3118 HypreLOBPCG::HypreMultiVector::GetVector(unsigned int i)
3119 {
3120  MFEM_ASSERT((int)i < nv, "index out of range");
3121 
3122  return ( *hpv[i] );
3123 }
3124 
3125 HypreParVector **
3126 HypreLOBPCG::HypreMultiVector::StealVectors()
3127 {
3128  HypreParVector ** hpv_ret = hpv;
3129 
3130  hpv = NULL;
3131 
3132  mv_TempMultiVector * mv_tmp =
3133  (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
3134 
3135  mv_tmp->ownsVectors = 0;
3136 
3137  for (int i=0; i<nv; i++)
3138  {
3139  hpv_ret[i]->SetOwnership(1);
3140  }
3141 
3142  return hpv_ret;
3143 }
3144 
3146  : comm(c),
3147  myid(0),
3148  numProcs(1),
3149  nev(10),
3150  seed(75),
3151  glbSize(-1),
3152  part(NULL),
3153  multi_vec(NULL),
3154  x(NULL),
3155  subSpaceProj(NULL)
3156 {
3157  MPI_Comm_size(comm,&numProcs);
3158  MPI_Comm_rank(comm,&myid);
3159 
3160  HYPRE_ParCSRSetupInterpreter(&interpreter);
3161  HYPRE_ParCSRSetupMatvec(&matvec_fn);
3162  HYPRE_LOBPCGCreate(&interpreter, &matvec_fn, &lobpcg_solver);
3163 }
3164 
3166 {
3167  delete multi_vec;
3168  delete x;
3169  delete [] part;
3170 
3171  HYPRE_LOBPCGDestroy(lobpcg_solver);
3172 }
3173 
3174 void
3176 {
3177  HYPRE_LOBPCGSetTol(lobpcg_solver, tol);
3178 }
3179 
3180 void
3181 HypreLOBPCG::SetRelTol(double rel_tol)
3182 {
3183 #if MFEM_HYPRE_VERSION >= 21101
3184  HYPRE_LOBPCGSetRTol(lobpcg_solver, rel_tol);
3185 #else
3186  MFEM_ABORT("This method requires HYPRE version >= 2.11.1");
3187 #endif
3188 }
3189 
3190 void
3192 {
3193  HYPRE_LOBPCGSetMaxIter(lobpcg_solver, max_iter);
3194 }
3195 
3196 void
3198 {
3199  if (myid == 0)
3200  {
3201  HYPRE_LOBPCGSetPrintLevel(lobpcg_solver, logging);
3202  }
3203 }
3204 
3205 void
3207 {
3208  HYPRE_LOBPCGSetPrecondUsageMode(lobpcg_solver, pcg_mode);
3209 }
3210 
3211 void
3213 {
3214  HYPRE_LOBPCGSetPrecond(lobpcg_solver,
3215  (HYPRE_PtrToSolverFcn)this->PrecondSolve,
3216  (HYPRE_PtrToSolverFcn)this->PrecondSetup,
3217  (HYPRE_Solver)&precond);
3218 }
3219 
3220 void
3222 {
3223  int locSize = A.Width();
3224 
3225  if (HYPRE_AssumedPartitionCheck())
3226  {
3227  part = new HYPRE_Int[2];
3228 
3229  MPI_Scan(&locSize, &part[1], 1, HYPRE_MPI_INT, MPI_SUM, comm);
3230 
3231  part[0] = part[1] - locSize;
3232 
3233  MPI_Allreduce(&locSize, &glbSize, 1, HYPRE_MPI_INT, MPI_SUM, comm);
3234  }
3235  else
3236  {
3237  part = new HYPRE_Int[numProcs+1];
3238 
3239  MPI_Allgather(&locSize, 1, MPI_INT,
3240  &part[1], 1, HYPRE_MPI_INT, comm);
3241 
3242  part[0] = 0;
3243  for (int i=0; i<numProcs; i++)
3244  {
3245  part[i+1] += part[i];
3246  }
3247 
3248  glbSize = part[numProcs];
3249  }
3250 
3251  if ( x != NULL )
3252  {
3253  delete x;
3254  }
3255 
3256  // Create a distributed vector without a data array.
3257  x = new HypreParVector(comm,glbSize,NULL,part);
3258 
3259  matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3260  matvec_fn.Matvec = this->OperatorMatvec;
3261  matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3262 
3263  HYPRE_LOBPCGSetup(lobpcg_solver,(HYPRE_Matrix)&A,NULL,NULL);
3264 }
3265 
3266 void
3268 {
3269  matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3270  matvec_fn.Matvec = this->OperatorMatvec;
3271  matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3272 
3273  HYPRE_LOBPCGSetupB(lobpcg_solver,(HYPRE_Matrix)&M,NULL);
3274 }
3275 
3276 void
3278 {
3279  // Initialize eigenvalues array with marker values
3280  eigs.SetSize(nev);
3281 
3282  for (int i=0; i<nev; i++)
3283  {
3284  eigs[i] = eigenvalues[i];
3285  }
3286 }
3287 
3290 {
3291  return multi_vec->GetVector(i);
3292 }
3293 
3294 void
3296 {
3297  // Initialize HypreMultiVector object if necessary
3298  if ( multi_vec == NULL )
3299  {
3300  MFEM_ASSERT(x != NULL, "In HypreLOBPCG::SetInitialVectors()");
3301 
3302  multi_vec = new HypreMultiVector(nev, *x, interpreter);
3303  }
3304 
3305  // Copy the vectors provided
3306  for (int i=0; i < min(num_vecs,nev); i++)
3307  {
3308  multi_vec->GetVector(i) = *vecs[i];
3309  }
3310 
3311  // Randomize any remaining vectors
3312  for (int i=min(num_vecs,nev); i < nev; i++)
3313  {
3314  multi_vec->GetVector(i).Randomize(seed);
3315  }
3316 
3317  // Ensure all vectors are in the proper subspace
3318  if ( subSpaceProj != NULL )
3319  {
3320  HypreParVector y(*x);
3321  y = multi_vec->GetVector(0);
3322 
3323  for (int i=1; i<nev; i++)
3324  {
3325  subSpaceProj->Mult(multi_vec->GetVector(i),
3326  multi_vec->GetVector(i-1));
3327  }
3328  subSpaceProj->Mult(y,
3329  multi_vec->GetVector(nev-1));
3330  }
3331 }
3332 
3333 void
3335 {
3336  // Initialize HypreMultiVector object if necessary
3337  if ( multi_vec == NULL )
3338  {
3339  MFEM_ASSERT(x != NULL, "In HypreLOBPCG::Solve()");
3340 
3341  multi_vec = new HypreMultiVector(nev, *x, interpreter);
3342  multi_vec->Randomize(seed);
3343 
3344  if ( subSpaceProj != NULL )
3345  {
3346  HypreParVector y(*x);
3347  y = multi_vec->GetVector(0);
3348 
3349  for (int i=1; i<nev; i++)
3350  {
3351  subSpaceProj->Mult(multi_vec->GetVector(i),
3352  multi_vec->GetVector(i-1));
3353  }
3354  subSpaceProj->Mult(y, multi_vec->GetVector(nev-1));
3355  }
3356  }
3357 
3358  eigenvalues.SetSize(nev);
3359  eigenvalues = NAN;
3360 
3361  // Perform eigenmode calculation
3362  //
3363  // The eigenvalues are computed in ascending order (internally the
3364  // order is determined by the LAPACK routine 'dsydv'.)
3365  HYPRE_LOBPCGSolve(lobpcg_solver, NULL, *multi_vec, eigenvalues);
3366 }
3367 
3368 void *
3369 HypreLOBPCG::OperatorMatvecCreate( void *A,
3370  void *x )
3371 {
3372  void *matvec_data;
3373 
3374  matvec_data = NULL;
3375 
3376  return ( matvec_data );
3377 }
3378 
3379 HYPRE_Int
3380 HypreLOBPCG::OperatorMatvec( void *matvec_data,
3381  HYPRE_Complex alpha,
3382  void *A,
3383  void *x,
3384  HYPRE_Complex beta,
3385  void *y )
3386 {
3387  MFEM_VERIFY(alpha == 1.0 && beta == 0.0, "values not supported");
3388 
3389  Operator *Aop = (Operator*)A;
3390 
3391  int width = Aop->Width();
3392 
3393  hypre_ParVector * xPar = (hypre_ParVector *)x;
3394  hypre_ParVector * yPar = (hypre_ParVector *)y;
3395 
3396  Vector xVec(xPar->local_vector->data, width);
3397  Vector yVec(yPar->local_vector->data, width);
3398 
3399  Aop->Mult( xVec, yVec );
3400 
3401  return 0;
3402 }
3403 
3404 HYPRE_Int
3405 HypreLOBPCG::OperatorMatvecDestroy( void *matvec_data )
3406 {
3407  return 0;
3408 }
3409 
3410 HYPRE_Int
3411 HypreLOBPCG::PrecondSolve(void *solver,
3412  void *A,
3413  void *b,
3414  void *x)
3415 {
3416  Solver *PC = (Solver*)solver;
3417  Operator *OP = (Operator*)A;
3418 
3419  int width = OP->Width();
3420 
3421  hypre_ParVector * bPar = (hypre_ParVector *)b;
3422  hypre_ParVector * xPar = (hypre_ParVector *)x;
3423 
3424  Vector bVec(bPar->local_vector->data, width);
3425  Vector xVec(xPar->local_vector->data, width);
3426 
3427  PC->Mult( bVec, xVec );
3428 
3429  return 0;
3430 }
3431 
3432 HYPRE_Int
3433 HypreLOBPCG::PrecondSetup(void *solver,
3434  void *A,
3435  void *b,
3436  void *x)
3437 {
3438  return 0;
3439 }
3440 
3441 HypreAME::HypreAME(MPI_Comm comm)
3442  : myid(0),
3443  numProcs(1),
3444  nev(10),
3445  setT(false),
3446  ams_precond(NULL),
3447  eigenvalues(NULL),
3448  multi_vec(NULL),
3449  eigenvectors(NULL)
3450 {
3451  MPI_Comm_size(comm,&numProcs);
3452  MPI_Comm_rank(comm,&myid);
3453 
3454  HYPRE_AMECreate(&ame_solver);
3455  HYPRE_AMESetPrintLevel(ame_solver, 0);
3456 }
3457 
3459 {
3460  if ( multi_vec )
3461  {
3462  hypre_TFree(multi_vec);
3463  }
3464 
3465  if ( eigenvectors )
3466  {
3467  for (int i=0; i<nev; i++)
3468  {
3469  delete eigenvectors[i];
3470  }
3471  }
3472  delete [] eigenvectors;
3473 
3474  if ( eigenvalues )
3475  {
3476  hypre_TFree(eigenvalues);
3477  }
3478 
3479  HYPRE_AMEDestroy(ame_solver);
3480 }
3481 
3482 void
3484 {
3485  nev = num_eigs;
3486 
3487  HYPRE_AMESetBlockSize(ame_solver, nev);
3488 }
3489 
3490 void
3491 HypreAME::SetTol(double tol)
3492 {
3493  HYPRE_AMESetTol(ame_solver, tol);
3494 }
3495 
3496 void
3497 HypreAME::SetRelTol(double rel_tol)
3498 {
3499 #if MFEM_HYPRE_VERSION >= 21101
3500  HYPRE_AMESetRTol(ame_solver, rel_tol);
3501 #else
3502  MFEM_ABORT("This method requires HYPRE version >= 2.11.1");
3503 #endif
3504 }
3505 
3506 void
3508 {
3509  HYPRE_AMESetMaxIter(ame_solver, max_iter);
3510 }
3511 
3512 void
3514 {
3515  if (myid == 0)
3516  {
3517  HYPRE_AMESetPrintLevel(ame_solver, logging);
3518  }
3519 }
3520 
3521 void
3523 {
3524  ams_precond = &precond;
3525 }
3526 
3527 void
3529 {
3530  if ( !setT )
3531  {
3532  HYPRE_Solver ams_precond_ptr = (HYPRE_Solver)*ams_precond;
3533 
3534  ams_precond->SetupFcn()(*ams_precond,A,NULL,NULL);
3535 
3536  HYPRE_AMESetAMSSolver(ame_solver, ams_precond_ptr);
3537  }
3538 
3539  HYPRE_AMESetup(ame_solver);
3540 }
3541 
3542 void
3544 {
3545  HYPRE_ParCSRMatrix parcsr_M = M;
3546  HYPRE_AMESetMassMatrix(ame_solver,(HYPRE_ParCSRMatrix)parcsr_M);
3547 }
3548 
3549 void
3551 {
3552  HYPRE_AMESolve(ame_solver);
3553 }
3554 
3555 void
3557 {
3558  // Initialize eigenvalues array with marker values
3559  eigs.SetSize(nev); eigs = -1.0;
3560 
3561  if ( eigenvalues == NULL )
3562  {
3563  // Grab eigenvalues from AME
3564  HYPRE_AMEGetEigenvalues(ame_solver,&eigenvalues);
3565  }
3566 
3567  // Copy eigenvalues to eigs array
3568  for (int i=0; i<nev; i++)
3569  {
3570  eigs[i] = eigenvalues[i];
3571  }
3572 }
3573 
3574 void
3575 HypreAME::createDummyVectors()
3576 {
3577  if ( multi_vec == NULL )
3578  {
3579  HYPRE_AMEGetEigenvectors(ame_solver,&multi_vec);
3580  }
3581 
3582  eigenvectors = new HypreParVector*[nev];
3583  for (int i=0; i<nev; i++)
3584  {
3585  eigenvectors[i] = new HypreParVector(multi_vec[i]);
3586  eigenvectors[i]->SetOwnership(1);
3587  }
3588 
3589 }
3590 
3591 HypreParVector &
3593 {
3594  if ( eigenvectors == NULL )
3595  {
3596  this->createDummyVectors();
3597  }
3598 
3599  return *eigenvectors[i];
3600 }
3601 
3602 HypreParVector **
3604 {
3605  if ( eigenvectors == NULL )
3606  {
3607  this->createDummyVectors();
3608  }
3609 
3610  // Set the local pointers to NULL so that they won't be deleted later
3611  HypreParVector ** vecs = eigenvectors;
3612  eigenvectors = NULL;
3613  multi_vec = NULL;
3614 
3615  return vecs;
3616 }
3617 
3618 }
3619 
3620 #endif
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:2620
void SetTol(double tol)
Definition: hypre.cpp:2088
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2251
int Size() const
Logical size of the array.
Definition: array.hpp:110
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1286
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:972
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
Definition: vector.hpp:51
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:298
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:127
int * GetJ()
Definition: table.hpp:111
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:514
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3289
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:544
int NumCols() const
Definition: array.hpp:287
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3070
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
Definition: hypre.cpp:1309
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1476
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:617
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:770
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
Definition: hypre.cpp:1329
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Definition: hypre.cpp:1023
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1017
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:320
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
Definition: hypre.hpp:421
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:614
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:546
void MakeDataOwner()
Definition: vector.hpp:107
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:42
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
Definition: hypre.cpp:3295
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:2830
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:1760
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:62
T * GetData()
Returns the data.
Definition: array.hpp:92
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3556
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:3522
int Size() const
Returns the size of the vector.
Definition: vector.hpp:113
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARHYP.
Definition: hypre.hpp:156
void SetTol(double tol)
Definition: hypre.cpp:3491
virtual ~HypreGMRES()
Definition: hypre.cpp:2327
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:3212
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:3267
Abstract parallel finite element space.
Definition: pfespace.hpp:31
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:707
void SetPrintLevel(int logging)
Definition: hypre.cpp:3197
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:263
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual ~HypreAMS()
Definition: hypre.cpp:2810
void SetTol(double tol)
Definition: hypre.cpp:2226
void SetOperator(HypreParMatrix &A)
Definition: hypre.cpp:3528
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition: sparsemat.hpp:401
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
Definition: hypre.hpp:540
double * GetData() const
Definition: vector.hpp:121
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:947
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2246
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1112
int Size_of_connections() const
Definition: table.hpp:95
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:528
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A * B.
Definition: hypre.cpp:1456
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2103
void CopyRowStarts()
Definition: hypre.cpp:794
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:530
void AddDomainInterpolator(DiscreteInterpolator *di)
void SetTol(double tol)
Definition: hypre.cpp:3175
HypreParVector * X1
Definition: hypre.hpp:514
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:1775
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:2259
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:169
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1530
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2502
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:3145
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2942
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:1899
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:184
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:374
int dim
Definition: ex3.cpp:47
void SetSymmetry(int sym)
Definition: hypre.cpp:2356
virtual ~HypreParaSails()
Definition: hypre.cpp:2361
MPI_Comm GetComm() const
Definition: pfespace.hpp:166
HYPRE_Int GetGlobalNumCols() const
Definition: hypre.hpp:377
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3191
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:538
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
void SetLogging(int logging)
Definition: hypre.cpp:2098
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:136
virtual ~HypreSolver()
Definition: hypre.cpp:2070
HyprePCG(HypreParMatrix &_A)
Definition: hypre.cpp:2077
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3497
void SetKDim(int dim)
Definition: hypre.cpp:2236
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
Definition: hypre.cpp:2116
int to_int(const std::string &str)
Definition: text.hpp:68
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:914
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:2023
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:178
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:3603
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:936
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:164
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:522
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2584
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:872
void Sort()
Sorts the array. This requires operator&lt; to be defined for T.
Definition: array.hpp:208
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:89
int Dimension() const
Definition: mesh.hpp:641
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:508
double * GetData() const
Return element data, i.e. array A.
Definition: sparsemat.hpp:135
HypreAMS(HypreParMatrix &A, ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:2631
int * GetI() const
Return the array I.
Definition: sparsemat.hpp:131
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:1782
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:2241
virtual ~HypreADS()
Definition: hypre.cpp:3048
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:255
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2093
int ParCSRRelax_FIR(hypre_ParCSRMatrix *A, hypre_ParVector *f, double max_eig, int poly_order, double *fir_coeffs, hypre_ParVector *u, hypre_ParVector *x0, hypre_ParVector *x1, hypre_ParVector *x2, hypre_ParVector *x3)
Definition: hypre.cpp:1617
HypreParMatrix * ParallelAssemble() const
Returns the matrix &quot;assembled&quot; on the true dofs.
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:1740
int SpaceDimension() const
Definition: mesh.hpp:642
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
HypreParVector * X
Definition: hypre.hpp:614
double Norml1() const
Returns the l_1 norm of the vector.
Definition: vector.cpp:717
int ParCSRRelax_Taubin(hypre_ParCSRMatrix *A, hypre_ParVector *f, double lambda, double mu, int N, double max_eig, hypre_ParVector *u, hypre_ParVector *r)
Definition: hypre.cpp:1580
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3277
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2825
HypreParVector * X
Definition: hypre.hpp:510
void mfem_error(const char *msg)
Definition: error.cpp:107
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:503
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3550
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:56
void SetNumModes(int num_eigs)
Definition: hypre.cpp:3483
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:361
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:1467
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:542
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3181
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:3441
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3507
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:1859
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:611
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:549
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:1226
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:203
virtual ~HypreSmoother()
Definition: hypre.cpp:1989
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
Definition: hypre.cpp:49
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2485
void SetData(double *_data)
Definition: hypre.cpp:159
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1147
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:1752
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:510
void SetMassMatrix(HypreParMatrix &M)
Definition: hypre.cpp:3543
virtual ~HyprePCG()
Definition: hypre.cpp:2202
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:368
void SetOperator(Operator &A)
Definition: hypre.cpp:3221
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2108
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:1746
void CopyColStarts()
Definition: hypre.cpp:831
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:520
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:1188
void SetPrintLevel(int logging)
Definition: hypre.cpp:3513
double ParNormlp(const Vector &vec, double p, MPI_Comm comm)
Compute the l_p norm of the Vector which is split without overlap across the given communicator...
Definition: hypre.cpp:214
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:124
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Definition: sparsemat.hpp:395
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:607
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:137
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:175
const double alpha
Definition: ex15.cpp:337
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:297
Vector data type.
Definition: vector.hpp:41
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:177
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
Definition: hypre.cpp:1350
HypreGMRES(HypreParMatrix &_A)
Definition: hypre.cpp:2208
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
Definition: hypre.cpp:780
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:3206
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:611
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2231
void GetOffd(SparseMatrix &offd, HYPRE_Int *&cmap) const
Get the local off-diagonal block. NOTE: &#39;offd&#39; will not own any data.
Definition: hypre.cpp:908
HYPRE_Int * GetColStarts() const
Definition: hypre.hpp:382
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3592
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:524
HYPRE_Int * GetRowStarts() const
Definition: hypre.hpp:380
Base class for solvers.
Definition: operator.hpp:259
Class for parallel grid function.
Definition: pgridfunc.hpp:32
double * data
Definition: vector.hpp:46
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:64
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:3132
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:174
virtual void Assemble(int skip_zeros=1)
int * GetJ() const
Return the array J.
Definition: sparsemat.hpp:133
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2129
Class for parallel meshes.
Definition: pmesh.hpp:29
HypreParVector * Z
Definition: hypre.hpp:512
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3334
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
Definition: hypre.cpp:1314
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:512
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:1734
HypreParaSails(HypreParMatrix &A)
Definition: hypre.cpp:2333
int NumRows() const
Definition: array.hpp:286
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
int poly_order
Order of the smoothing polynomial.
Definition: hypre.hpp:526
double lambda
Taubin&#39;s lambda-mu method parameters.
Definition: hypre.hpp:533
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:152