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