MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
complex_operator.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "complex_operator.hpp"
13 #include <set>
14 #include <map>
15 
16 // Define macro wrappers for hypre_TAlloc, hypre_CTAlloc and hypre_TFree:
17 // mfem_hypre_TAlloc, mfem_hypre_CTAlloc, and mfem_hypre_TFree, respectively.
18 // Note: the same macros are defined in hypre.cpp and hypre_parser.cpp.
19 #if MFEM_HYPRE_VERSION < 21400
20 
21 #define mfem_hypre_TAlloc(type, size) hypre_TAlloc(type, size)
22 #define mfem_hypre_CTAlloc(type, size) hypre_CTAlloc(type, size)
23 #define mfem_hypre_TFree(ptr) hypre_TFree(ptr)
24 
25 #else // MFEM_HYPRE_VERSION >= 21400
26 
27 // See the notes about hypre 2.14.0 in hypre.cpp
28 #define mfem_hypre_TAlloc(type, size) \
29  hypre_TAlloc(type, size, HYPRE_MEMORY_HOST)
30 #define mfem_hypre_CTAlloc(type, size) \
31  hypre_CTAlloc(type, size, HYPRE_MEMORY_HOST)
32 #define mfem_hypre_TFree(ptr) hypre_TFree(ptr, HYPRE_MEMORY_HOST)
33 
34 #endif // #if MFEM_HYPRE_VERSION < 21400
35 
36 namespace mfem
37 {
38 
40  bool ownReal, bool ownImag,
41  Convention convention)
42  : Operator(2*((Op_Real)?Op_Real->Height():Op_Imag->Height()),
43  2*((Op_Real)?Op_Real->Width():Op_Imag->Width()))
44  , Op_Real_(Op_Real)
45  , Op_Imag_(Op_Imag)
46  , ownReal_(ownReal)
47  , ownImag_(ownImag)
48  , convention_(convention)
49  , x_r_(NULL, width / 2)
50  , x_i_(NULL, width / 2)
51  , y_r_(NULL, height / 2)
52  , y_i_(NULL, height / 2)
53  , u_(NULL)
54  , v_(NULL)
55 {}
56 
58 {
59  if (ownReal_) { delete Op_Real_; }
60  if (ownImag_) { delete Op_Imag_; }
61  delete u_;
62  delete v_;
63 }
64 
66 {
67  MFEM_ASSERT(Op_Real_, "ComplexOperator has no real part!");
68  return *Op_Real_;
69 }
70 
72 {
73  MFEM_ASSERT(Op_Imag_, "ComplexOperator has no imaginary part!");
74  return *Op_Imag_;
75 }
76 
78 {
79  MFEM_ASSERT(Op_Real_, "ComplexOperator has no real part!");
80  return *Op_Real_;
81 }
82 
84 {
85  MFEM_ASSERT(Op_Imag_, "ComplexOperator has no imaginary part!");
86  return *Op_Imag_;
87 }
88 
89 void ComplexOperator::Mult(const Vector &x, Vector &y) const
90 {
91  double * x_data = x.GetData();
92  x_r_.SetData(x_data);
93  x_i_.SetData(&x_data[width / 2]);
94 
95  y_r_.SetData(&y[0]);
96  y_i_.SetData(&y[height / 2]);
97 
98  this->Mult(x_r_, x_i_, y_r_, y_i_);
99 }
100 
101 void ComplexOperator::Mult(const Vector &x_r, const Vector &x_i,
102  Vector &y_r, Vector &y_i) const
103 {
104  if (Op_Real_)
105  {
106  Op_Real_->Mult(x_r, y_r);
107  Op_Real_->Mult(x_i, y_i);
108  }
109  else
110  {
111  y_r = 0.0;
112  y_i = 0.0;
113  }
114  if (Op_Imag_)
115  {
116  if (!v_) { v_ = new Vector(Op_Imag_->Height()); }
117  Op_Imag_->Mult(x_i, *v_);
118  y_r_ -= *v_;
119  Op_Imag_->Mult(x_r, *v_);
120  y_i_ += *v_;
121  }
122 
124  {
125  y_i_ *= -1.0;
126  }
127 }
128 
130 {
131  double * x_data = x.GetData();
132  y_r_.SetData(x_data);
133  y_i_.SetData(&x_data[height / 2]);
134 
135  x_r_.SetData(&y[0]);
136  x_i_.SetData(&y[width / 2]);
137 
138  this->MultTranspose(y_r_, y_i_, x_r_, x_i_);
139 }
140 
141 void ComplexOperator::MultTranspose(const Vector &x_r, const Vector &x_i,
142  Vector &y_r, Vector &y_i) const
143 {
144  if (Op_Real_)
145  {
146  Op_Real_->MultTranspose(x_r, y_r);
147  Op_Real_->MultTranspose(x_i, y_i);
148 
150  {
151  y_i *= -1.0;
152  }
153  }
154  else
155  {
156  y_r = 0.0;
157  y_i = 0.0;
158  }
159  if (Op_Imag_)
160  {
161  if (!u_) { u_ = new Vector(Op_Imag_->Width()); }
162  Op_Imag_->MultTranspose(x_i, *u_);
163  y_r_.Add(convention_ == BLOCK_SYMMETRIC ? -1.0 : 1.0, *u_);
164  Op_Imag_->MultTranspose(x_r, *u_);
165  y_i_ -= *u_;
166  }
167 }
168 
169 
171 {
172  MFEM_ASSERT(Op_Real_, "ComplexSparseMatrix has no real part!");
173  return dynamic_cast<SparseMatrix &>(*Op_Real_);
174 }
175 
177 {
178  MFEM_ASSERT(Op_Imag_, "ComplexSparseMatrix has no imaginary part!");
179  return dynamic_cast<SparseMatrix &>(*Op_Imag_);
180 }
181 
183 {
184  MFEM_ASSERT(Op_Real_, "ComplexSparseMatrix has no real part!");
185  return dynamic_cast<const SparseMatrix &>(*Op_Real_);
186 }
187 
189 {
190  MFEM_ASSERT(Op_Imag_, "ComplexSparseMatrix has no imaginary part!");
191  return dynamic_cast<const SparseMatrix &>(*Op_Imag_);
192 }
193 
195 {
196  SparseMatrix * A_r = dynamic_cast<SparseMatrix*>(Op_Real_);
197  SparseMatrix * A_i = dynamic_cast<SparseMatrix*>(Op_Imag_);
198 
199  const int nrows_r = (A_r)?A_r->Height():0;
200  const int nrows_i = (A_i)?A_i->Height():0;
201  const int nrows = std::max(nrows_r, nrows_i);
202 
203  const int *I_r = (A_r)?A_r->GetI():NULL;
204  const int *I_i = (A_i)?A_i->GetI():NULL;
205 
206  const int *J_r = (A_r)?A_r->GetJ():NULL;
207  const int *J_i = (A_i)?A_i->GetJ():NULL;
208 
209  const double *D_r = (A_r)?A_r->GetData():NULL;
210  const double *D_i = (A_i)?A_i->GetData():NULL;
211 
212  const int nnz_r = (I_r)?I_r[nrows]:0;
213  const int nnz_i = (I_i)?I_i[nrows]:0;
214  const int nnz = 2 * (nnz_r + nnz_i);
215 
216  int *I = Memory<int>(this->Height()+1);
217  int *J = Memory<int>(nnz);
218  double *D = Memory<double>(nnz);
219 
220  const double factor = (convention_ == HERMITIAN) ? 1.0 : -1.0;
221 
222  I[0] = 0;
223  I[nrows] = nnz_r + nnz_i;
224  for (int i=0; i<nrows; i++)
225  {
226  I[i + 1] = ((I_r)?I_r[i+1]:0) + ((I_i)?I_i[i+1]:0);
227  I[i + nrows + 1] = I[i+1] + nnz_r + nnz_i;
228 
229  if (I_r)
230  {
231  const int off_i = (I_i)?(I_i[i+1] - I_i[i]):0;
232  for (int j=0; j<I_r[i+1] - I_r[i]; j++)
233  {
234  J[I[i] + j] = J_r[I_r[i] + j];
235  D[I[i] + j] = D_r[I_r[i] + j];
236 
237  J[I[i+nrows] + off_i + j] = J_r[I_r[i] + j] + nrows;
238  D[I[i+nrows] + off_i + j] = factor*D_r[I_r[i] + j];
239  }
240  }
241  if (I_i)
242  {
243  const int off_r = (I_r)?(I_r[i+1] - I_r[i]):0;
244  for (int j=0; j<I_i[i+1] - I_i[i]; j++)
245  {
246  J[I[i] + off_r + j] = J_i[I_i[i] + j] + nrows;
247  D[I[i] + off_r + j] = -D_i[I_i[i] + j];
248 
249  J[I[i+nrows] + j] = J_i[I_i[i] + j];
250  D[I[i+nrows] + j] = factor*D_i[I_i[i] + j];
251  }
252  }
253  }
254 
255  return new SparseMatrix(I, J, D, this->Height(), this->Width());
256 }
257 
258 #ifdef MFEM_USE_MPI
259 
261  HypreParMatrix * A_Imag,
262  bool ownReal, bool ownImag,
263  Convention convention)
264  : ComplexOperator(A_Real, A_Imag, ownReal, ownImag, convention)
265 {
266  comm_ = (A_Real) ? A_Real->GetComm() :
267  ((A_Imag) ? A_Imag->GetComm() : MPI_COMM_WORLD);
268 
269  MPI_Comm_rank(comm_, &myid_);
270  MPI_Comm_size(comm_, &nranks_);
271 }
272 
274 {
275  MFEM_ASSERT(Op_Real_, "ComplexHypreParMatrix has no real part!");
276  return dynamic_cast<HypreParMatrix &>(*Op_Real_);
277 }
278 
280 {
281  MFEM_ASSERT(Op_Imag_, "ComplexHypreParMatrix has no imaginary part!");
282  return dynamic_cast<HypreParMatrix &>(*Op_Imag_);
283 }
284 
286 {
287  MFEM_ASSERT(Op_Real_, "ComplexHypreParMatrix has no real part!");
288  return dynamic_cast<const HypreParMatrix &>(*Op_Real_);
289 }
290 
292 {
293  MFEM_ASSERT(Op_Imag_, "ComplexHypreParMatrix has no imaginary part!");
294  return dynamic_cast<const HypreParMatrix &>(*Op_Imag_);
295 }
296 
298 {
299  HypreParMatrix * A_r = dynamic_cast<HypreParMatrix*>(Op_Real_);
300  HypreParMatrix * A_i = dynamic_cast<HypreParMatrix*>(Op_Imag_);
301 
302  if ( A_r == NULL && A_i == NULL ) { return NULL; }
303 
304  HYPRE_Int global_num_rows_r = (A_r) ? A_r->GetGlobalNumRows() : 0;
305  HYPRE_Int global_num_rows_i = (A_i) ? A_i->GetGlobalNumRows() : 0;
306  HYPRE_Int global_num_rows = std::max(global_num_rows_r, global_num_rows_i);
307 
308  HYPRE_Int global_num_cols_r = (A_r) ? A_r->GetGlobalNumCols() : 0;
309  HYPRE_Int global_num_cols_i = (A_i) ? A_i->GetGlobalNumCols() : 0;
310  HYPRE_Int global_num_cols = std::max(global_num_cols_r, global_num_cols_i);
311 
312  int row_starts_size = (HYPRE_AssumedPartitionCheck()) ? 2 : nranks_ + 1;
313  HYPRE_Int * row_starts = mfem_hypre_CTAlloc(HYPRE_Int, row_starts_size);
314  HYPRE_Int * col_starts = mfem_hypre_CTAlloc(HYPRE_Int, row_starts_size);
315 
316  const HYPRE_Int * row_starts_z = (A_r) ? A_r->RowPart() :
317  ((A_i) ? A_i->RowPart() : NULL);
318  const HYPRE_Int * col_starts_z = (A_r) ? A_r->ColPart() :
319  ((A_i) ? A_i->ColPart() : NULL);
320 
321  for (int i = 0; i < row_starts_size; i++)
322  {
323  row_starts[i] = 2 * row_starts_z[i];
324  col_starts[i] = 2 * col_starts_z[i];
325  }
326 
327  SparseMatrix diag_r, diag_i, offd_r, offd_i;
328  HYPRE_Int * cmap_r, * cmap_i;
329 
330  int nrows_r = 0, nrows_i = 0, ncols_r = 0, ncols_i = 0;
331  int ncols_offd_r = 0, ncols_offd_i = 0;
332  if (A_r)
333  {
334  A_r->GetDiag(diag_r);
335  A_r->GetOffd(offd_r, cmap_r);
336  nrows_r = diag_r.Height();
337  ncols_r = diag_r.Width();
338  ncols_offd_r = offd_r.Width();
339  }
340  if (A_i)
341  {
342  A_i->GetDiag(diag_i);
343  A_i->GetOffd(offd_i, cmap_i);
344  nrows_i = diag_i.Height();
345  ncols_i = diag_i.Width();
346  ncols_offd_i = offd_i.Width();
347  }
348  int nrows = std::max(nrows_r, nrows_i);
349  int ncols = std::max(ncols_r, ncols_i);
350 
351  // Determine the unique set of off-diagonal columns global indices
352  std::set<int> cset;
353  for (int i=0; i<ncols_offd_r; i++)
354  {
355  cset.insert(cmap_r[i]);
356  }
357  for (int i=0; i<ncols_offd_i; i++)
358  {
359  cset.insert(cmap_i[i]);
360  }
361  int num_cols_offd = (int)cset.size();
362 
363  // Extract pointers to the various CSR arrays of the diagonal blocks
364  const int * diag_r_I = (A_r) ? diag_r.GetI() : NULL;
365  const int * diag_i_I = (A_i) ? diag_i.GetI() : NULL;
366 
367  const int * diag_r_J = (A_r) ? diag_r.GetJ() : NULL;
368  const int * diag_i_J = (A_i) ? diag_i.GetJ() : NULL;
369 
370  const double * diag_r_D = (A_r) ? diag_r.GetData() : NULL;
371  const double * diag_i_D = (A_i) ? diag_i.GetData() : NULL;
372 
373  int diag_r_nnz = (diag_r_I) ? diag_r_I[nrows] : 0;
374  int diag_i_nnz = (diag_i_I) ? diag_i_I[nrows] : 0;
375  int diag_nnz = 2 * (diag_r_nnz + diag_i_nnz);
376 
377  // Extract pointers to the various CSR arrays of the off-diagonal blocks
378  const int * offd_r_I = (A_r) ? offd_r.GetI() : NULL;
379  const int * offd_i_I = (A_i) ? offd_i.GetI() : NULL;
380 
381  const int * offd_r_J = (A_r) ? offd_r.GetJ() : NULL;
382  const int * offd_i_J = (A_i) ? offd_i.GetJ() : NULL;
383 
384  const double * offd_r_D = (A_r) ? offd_r.GetData() : NULL;
385  const double * offd_i_D = (A_i) ? offd_i.GetData() : NULL;
386 
387  int offd_r_nnz = (offd_r_I) ? offd_r_I[nrows] : 0;
388  int offd_i_nnz = (offd_i_I) ? offd_i_I[nrows] : 0;
389  int offd_nnz = 2 * (offd_r_nnz + offd_i_nnz);
390 
391  // Allocate CSR arrays for the combined matrix
392  HYPRE_Int * diag_I = mfem_hypre_CTAlloc(HYPRE_Int, 2 * nrows + 1);
393  HYPRE_Int * diag_J = mfem_hypre_CTAlloc(HYPRE_Int, diag_nnz);
394  double * diag_D = mfem_hypre_CTAlloc(double, diag_nnz);
395 
396  HYPRE_Int * offd_I = mfem_hypre_CTAlloc(HYPRE_Int, 2 * nrows + 1);
397  HYPRE_Int * offd_J = mfem_hypre_CTAlloc(HYPRE_Int, offd_nnz);
398  double * offd_D = mfem_hypre_CTAlloc(double, offd_nnz);
399  HYPRE_Int * cmap = mfem_hypre_CTAlloc(HYPRE_Int, 2 * num_cols_offd);
400 
401  // Fill the CSR arrays for the diagonal portion of the matrix
402  const double factor = (convention_ == HERMITIAN) ? 1.0 : -1.0;
403 
404  diag_I[0] = 0;
405  diag_I[nrows] = diag_r_nnz + diag_i_nnz;
406  for (int i=0; i<nrows; i++)
407  {
408  diag_I[i + 1] = ((diag_r_I)?diag_r_I[i+1]:0) +
409  ((diag_i_I)?diag_i_I[i+1]:0);
410  diag_I[i + nrows + 1] = diag_I[i+1] + diag_r_nnz + diag_i_nnz;
411 
412  if (diag_r_I)
413  {
414  for (int j=0; j<diag_r_I[i+1] - diag_r_I[i]; j++)
415  {
416  diag_J[diag_I[i] + j] = diag_r_J[diag_r_I[i] + j];
417  diag_D[diag_I[i] + j] = diag_r_D[diag_r_I[i] + j];
418 
419  diag_J[diag_I[i+nrows] + j] =
420  diag_r_J[diag_r_I[i] + j] + ncols;
421  diag_D[diag_I[i+nrows] + j] =
422  factor * diag_r_D[diag_r_I[i] + j];
423  }
424  }
425  if (diag_i_I)
426  {
427  const int off_r = (diag_r_I)?(diag_r_I[i+1] - diag_r_I[i]):0;
428  for (int j=0; j<diag_i_I[i+1] - diag_i_I[i]; j++)
429  {
430  diag_J[diag_I[i] + off_r + j] = diag_i_J[diag_i_I[i] + j] + ncols;
431  diag_D[diag_I[i] + off_r + j] = -diag_i_D[diag_i_I[i] + j];
432 
433  diag_J[diag_I[i+nrows] + off_r + j] = diag_i_J[diag_i_I[i] + j];
434  diag_D[diag_I[i+nrows] + off_r + j] =
435  factor * diag_i_D[diag_i_I[i] + j];
436  }
437  }
438  }
439 
440  // Determine the mappings describing the layout of off-diagonal columns
441  int num_recv_procs = 0;
442  HYPRE_Int * offd_col_start_stop = NULL;
443  this->getColStartStop(A_r, A_i, num_recv_procs, offd_col_start_stop);
444 
445  std::set<int>::iterator sit;
446  std::map<int,int> cmapa, cmapb, cinvmap;
447  for (sit=cset.begin(); sit!=cset.end(); sit++)
448  {
449  int col_orig = *sit;
450  int col_2x2 = -1;
451  int col_size = 0;
452  for (int i=0; i<num_recv_procs; i++)
453  {
454  if (offd_col_start_stop[2*i] <= col_orig &&
455  col_orig < offd_col_start_stop[2*i+1])
456  {
457  col_2x2 = offd_col_start_stop[2*i] + col_orig;
458  col_size = offd_col_start_stop[2*i+1] - offd_col_start_stop[2*i];
459  break;
460  }
461  }
462  cmapa[*sit] = col_2x2;
463  cmapb[*sit] = col_2x2 + col_size;
464  cinvmap[col_2x2] = -1;
465  cinvmap[col_2x2 + col_size] = -1;
466  }
467  delete [] offd_col_start_stop;
468 
469  std::map<int, int>::iterator mit;
470  int i = 0;
471  for (mit=cinvmap.begin(); mit!=cinvmap.end(); mit++, i++)
472  {
473  mit->second = i;
474  cmap[i] = mit->first;
475  }
476 
477  // Fill the CSR arrays for the off-diagonal portion of the matrix
478  offd_I[0] = 0;
479  offd_I[nrows] = offd_r_nnz + offd_i_nnz;
480  for (int i=0; i<nrows; i++)
481  {
482  offd_I[i + 1] = ((offd_r_I)?offd_r_I[i+1]:0) +
483  ((offd_i_I)?offd_i_I[i+1]:0);
484  offd_I[i + nrows + 1] = offd_I[i+1] + offd_r_nnz + offd_i_nnz;
485 
486  if (offd_r_I)
487  {
488  const int off_i = (offd_i_I)?(offd_i_I[i+1] - offd_i_I[i]):0;
489  for (int j=0; j<offd_r_I[i+1] - offd_r_I[i]; j++)
490  {
491  offd_J[offd_I[i] + j] =
492  cinvmap[cmapa[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
493  offd_D[offd_I[i] + j] = offd_r_D[offd_r_I[i] + j];
494 
495  offd_J[offd_I[i+nrows] + off_i + j] =
496  cinvmap[cmapb[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
497  offd_D[offd_I[i+nrows] + off_i + j] =
498  factor * offd_r_D[offd_r_I[i] + j];
499  }
500  }
501  if (offd_i_I)
502  {
503  const int off_r = (offd_r_I)?(offd_r_I[i+1] - offd_r_I[i]):0;
504  for (int j=0; j<offd_i_I[i+1] - offd_i_I[i]; j++)
505  {
506  offd_J[offd_I[i] + off_r + j] =
507  cinvmap[cmapb[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
508  offd_D[offd_I[i] + off_r + j] = -offd_i_D[offd_i_I[i] + j];
509 
510  offd_J[offd_I[i+nrows] + j] =
511  cinvmap[cmapa[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
512  offd_D[offd_I[i+nrows] + j] = factor * offd_i_D[offd_i_I[i] + j];
513  }
514  }
515  }
516 
517  // Construct the combined matrix
518  HypreParMatrix * A = new HypreParMatrix(comm_,
519  2 * global_num_rows,
520  2 * global_num_cols,
521  row_starts, col_starts,
522  diag_I, diag_J, diag_D,
523  offd_I, offd_J, offd_D,
524  2 * num_cols_offd, cmap);
525 
526  // Give the new matrix ownership of its internal arrays
527  A->SetOwnerFlags(-1,-1,-1);
528  hypre_CSRMatrixSetDataOwner(((hypre_ParCSRMatrix*)(*A))->diag,1);
529  hypre_CSRMatrixSetDataOwner(((hypre_ParCSRMatrix*)(*A))->offd,1);
530  hypre_ParCSRMatrixSetRowStartsOwner((hypre_ParCSRMatrix*)(*A),1);
531  hypre_ParCSRMatrixSetColStartsOwner((hypre_ParCSRMatrix*)(*A),1);
532 
533  return A;
534 }
535 
536 void
537 ComplexHypreParMatrix::getColStartStop(const HypreParMatrix * A_r,
538  const HypreParMatrix * A_i,
539  int & num_recv_procs,
540  HYPRE_Int *& offd_col_start_stop) const
541 {
542  hypre_ParCSRCommPkg * comm_pkg_r =
543  (A_r) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_r)) : NULL;
544  hypre_ParCSRCommPkg * comm_pkg_i =
545  (A_i) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_i)) : NULL;
546 
547  std::set<HYPRE_Int> send_procs, recv_procs;
548  if ( comm_pkg_r )
549  {
550  for (HYPRE_Int i=0; i<comm_pkg_r->num_sends; i++)
551  {
552  send_procs.insert(comm_pkg_r->send_procs[i]);
553  }
554  for (HYPRE_Int i=0; i<comm_pkg_r->num_recvs; i++)
555  {
556  recv_procs.insert(comm_pkg_r->recv_procs[i]);
557  }
558  }
559  if ( comm_pkg_i )
560  {
561  for (HYPRE_Int i=0; i<comm_pkg_i->num_sends; i++)
562  {
563  send_procs.insert(comm_pkg_i->send_procs[i]);
564  }
565  for (HYPRE_Int i=0; i<comm_pkg_i->num_recvs; i++)
566  {
567  recv_procs.insert(comm_pkg_i->recv_procs[i]);
568  }
569  }
570 
571  num_recv_procs = (int)recv_procs.size();
572 
573  HYPRE_Int loc_start_stop[2];
574  offd_col_start_stop = new HYPRE_Int[2 * num_recv_procs];
575 
576  const HYPRE_Int * row_part = (A_r) ? A_r->RowPart() :
577  ((A_i) ? A_i->RowPart() : NULL);
578 
579  int row_part_ind = (HYPRE_AssumedPartitionCheck()) ? 0 : myid_;
580  loc_start_stop[0] = row_part[row_part_ind];
581  loc_start_stop[1] = row_part[row_part_ind+1];
582 
583  MPI_Request * req = new MPI_Request[send_procs.size()+recv_procs.size()];
584  MPI_Status * stat = new MPI_Status[send_procs.size()+recv_procs.size()];
585  int send_count = 0;
586  int recv_count = 0;
587  int tag = 0;
588 
589  std::set<HYPRE_Int>::iterator sit;
590  for (sit=send_procs.begin(); sit!=send_procs.end(); sit++)
591  {
592  MPI_Isend(loc_start_stop, 2, HYPRE_MPI_INT,
593  *sit, tag, comm_, &req[send_count]);
594  send_count++;
595  }
596  for (sit=recv_procs.begin(); sit!=recv_procs.end(); sit++)
597  {
598  MPI_Irecv(&offd_col_start_stop[2*recv_count], 2, HYPRE_MPI_INT,
599  *sit, tag, comm_, &req[send_count+recv_count]);
600  recv_count++;
601  }
602 
603  MPI_Waitall(send_count+recv_count, req, stat);
604 
605  delete [] req;
606  delete [] stat;
607 }
608 
609 #endif // MFEM_USE_MPI
610 
611 }
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:336
virtual Operator & real()
Real or imaginary part accessor methods.
ComplexHypreParMatrix(HypreParMatrix *A_Real, HypreParMatrix *A_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:376
Mimic the action of a complex operator using two real operators.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:63
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:146
int * GetI()
Return the array I.
Definition: sparsemat.hpp:141
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Alternate convention for damping operators.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
ComplexOperator(Operator *Op_Real, Operator *Op_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
Constructs complex operator object.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:84
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:151
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:425
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:57
void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
Definition: hypre.hpp:348
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:939
void SetData(double *d)
Definition: vector.hpp:118
virtual SparseMatrix & real()
Real or imaginary part accessor methods.
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:372
SparseMatrix * GetSystemMatrix() const
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual Operator & imag()
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:190
virtual SparseMatrix & imag()
Native convention for Hermitian operators.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Vector data type.
Definition: vector.hpp:48
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:975
HypreParMatrix * GetSystemMatrix() const
virtual HypreParMatrix & imag()
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual HypreParMatrix & real()
Real or imaginary part accessor methods.
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28