MFEM  v4.2.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 namespace mfem
17 {
18 
20  bool ownReal, bool ownImag,
21  Convention convention)
22  : Operator(2*((Op_Real)?Op_Real->Height():Op_Imag->Height()),
23  2*((Op_Real)?Op_Real->Width():Op_Imag->Width()))
24  , Op_Real_(Op_Real)
25  , Op_Imag_(Op_Imag)
26  , ownReal_(ownReal)
27  , ownImag_(ownImag)
28  , convention_(convention)
29  , x_r_()
30  , x_i_()
31  , y_r_()
32  , y_i_()
33  , u_(NULL)
34  , v_(NULL)
35 {}
36 
38 {
39  if (ownReal_) { delete Op_Real_; }
40  if (ownImag_) { delete Op_Imag_; }
41  delete u_;
42  delete v_;
43 }
44 
46 {
47  MFEM_ASSERT(Op_Real_, "ComplexOperator has no real part!");
48  return *Op_Real_;
49 }
50 
52 {
53  MFEM_ASSERT(Op_Imag_, "ComplexOperator has no imaginary part!");
54  return *Op_Imag_;
55 }
56 
58 {
59  MFEM_ASSERT(Op_Real_, "ComplexOperator has no real part!");
60  return *Op_Real_;
61 }
62 
64 {
65  MFEM_ASSERT(Op_Imag_, "ComplexOperator has no imaginary part!");
66  return *Op_Imag_;
67 }
68 
69 void ComplexOperator::Mult(const Vector &x, Vector &y) const
70 {
71  x.Read();
72  y.UseDevice(true); y = 0.0;
73 
74  x_r_.MakeRef(const_cast<Vector&>(x), 0, width/2);
75  x_i_.MakeRef(const_cast<Vector&>(x), width/2, width/2);
76 
77  y_r_.MakeRef(y, 0, height/2);
78  y_i_.MakeRef(y, height/2, height/2);
79 
80  this->Mult(x_r_, x_i_, y_r_, y_i_);
81 
84 
85  // Destroy alias vectors to prevent dangling aliases when the base vectors
86  // are deleted
87  x_r_.Destroy();
88  x_i_.Destroy();
89  y_r_.Destroy();
90  y_i_.Destroy();
91 }
92 
93 void ComplexOperator::Mult(const Vector &x_r, const Vector &x_i,
94  Vector &y_r, Vector &y_i) const
95 {
96  if (Op_Real_)
97  {
98  Op_Real_->Mult(x_r, y_r);
99  Op_Real_->Mult(x_i, y_i);
100  }
101  else
102  {
103  y_r = 0.0;
104  y_i = 0.0;
105  }
106 
107  if (Op_Imag_)
108  {
109  if (!v_) { v_ = new Vector(); }
110  v_->UseDevice(true);
111  v_->SetSize(Op_Imag_->Height());
112 
113  Op_Imag_->Mult(x_i, *v_);
114  y_r.Add(-1.0, *v_);
115  Op_Imag_->Mult(x_r, *v_);
116  y_i.Add(1.0, *v_);
117  }
118 
120  {
121  y_i *= -1.0;
122  }
123 }
124 
126 {
127  x.Read();
128  y.UseDevice(true); y = 0.0;
129 
130  x_r_.MakeRef(const_cast<Vector&>(x), 0, height/2);
131  x_i_.MakeRef(const_cast<Vector&>(x), height/2, height/2);
132 
133  y_r_.MakeRef(y, 0, width/2);
134  y_i_.MakeRef(y, width/2, width/2);
135 
136  this->MultTranspose(x_r_, x_i_, y_r_, y_i_);
137 
140 
141  // Destroy alias vectors to prevent dangling aliases when the base vectors
142  // are deleted
143  x_r_.Destroy();
144  x_i_.Destroy();
145  y_r_.Destroy();
146  y_i_.Destroy();
147 }
148 
149 void ComplexOperator::MultTranspose(const Vector &x_r, const Vector &x_i,
150  Vector &y_r, Vector &y_i) const
151 {
152  if (Op_Real_)
153  {
154  Op_Real_->MultTranspose(x_r, y_r);
155  Op_Real_->MultTranspose(x_i, y_i);
156 
158  {
159  y_i *= -1.0;
160  }
161  }
162  else
163  {
164  y_r = 0.0;
165  y_i = 0.0;
166  }
167 
168  if (Op_Imag_)
169  {
170  if (!u_) { u_ = new Vector(); }
171  u_->UseDevice(true);
172  u_->SetSize(Op_Imag_->Width());
173 
174  Op_Imag_->MultTranspose(x_i, *u_);
175  y_r.Add(convention_ == BLOCK_SYMMETRIC ? -1.0 : 1.0, *u_);
176  Op_Imag_->MultTranspose(x_r, *u_);
177  y_i.Add(-1.0, *u_);
178  }
179 }
180 
181 
183 {
184  MFEM_ASSERT(Op_Real_, "ComplexSparseMatrix has no real part!");
185  return dynamic_cast<SparseMatrix &>(*Op_Real_);
186 }
187 
189 {
190  MFEM_ASSERT(Op_Imag_, "ComplexSparseMatrix has no imaginary part!");
191  return dynamic_cast<SparseMatrix &>(*Op_Imag_);
192 }
193 
195 {
196  MFEM_ASSERT(Op_Real_, "ComplexSparseMatrix has no real part!");
197  return dynamic_cast<const SparseMatrix &>(*Op_Real_);
198 }
199 
201 {
202  MFEM_ASSERT(Op_Imag_, "ComplexSparseMatrix has no imaginary part!");
203  return dynamic_cast<const SparseMatrix &>(*Op_Imag_);
204 }
205 
207 {
208  SparseMatrix * A_r = dynamic_cast<SparseMatrix*>(Op_Real_);
209  SparseMatrix * A_i = dynamic_cast<SparseMatrix*>(Op_Imag_);
210 
211  const int nrows_r = (A_r)?A_r->Height():0;
212  const int nrows_i = (A_i)?A_i->Height():0;
213  const int nrows = std::max(nrows_r, nrows_i);
214 
215  const int *I_r = (A_r)?A_r->GetI():NULL;
216  const int *I_i = (A_i)?A_i->GetI():NULL;
217 
218  const int *J_r = (A_r)?A_r->GetJ():NULL;
219  const int *J_i = (A_i)?A_i->GetJ():NULL;
220 
221  const double *D_r = (A_r)?A_r->GetData():NULL;
222  const double *D_i = (A_i)?A_i->GetData():NULL;
223 
224  const int nnz_r = (I_r)?I_r[nrows]:0;
225  const int nnz_i = (I_i)?I_i[nrows]:0;
226  const int nnz = 2 * (nnz_r + nnz_i);
227 
228  int *I = Memory<int>(this->Height()+1);
229  int *J = Memory<int>(nnz);
230  double *D = Memory<double>(nnz);
231 
232  const double factor = (convention_ == HERMITIAN) ? 1.0 : -1.0;
233 
234  I[0] = 0;
235  I[nrows] = nnz_r + nnz_i;
236  for (int i=0; i<nrows; i++)
237  {
238  I[i + 1] = ((I_r)?I_r[i+1]:0) + ((I_i)?I_i[i+1]:0);
239  I[i + nrows + 1] = I[i+1] + nnz_r + nnz_i;
240 
241  if (I_r)
242  {
243  const int off_i = (I_i)?(I_i[i+1] - I_i[i]):0;
244  for (int j=0; j<I_r[i+1] - I_r[i]; j++)
245  {
246  J[I[i] + j] = J_r[I_r[i] + j];
247  D[I[i] + j] = D_r[I_r[i] + j];
248 
249  J[I[i+nrows] + off_i + j] = J_r[I_r[i] + j] + nrows;
250  D[I[i+nrows] + off_i + j] = factor*D_r[I_r[i] + j];
251  }
252  }
253  if (I_i)
254  {
255  const int off_r = (I_r)?(I_r[i+1] - I_r[i]):0;
256  for (int j=0; j<I_i[i+1] - I_i[i]; j++)
257  {
258  J[I[i] + off_r + j] = J_i[I_i[i] + j] + nrows;
259  D[I[i] + off_r + j] = -D_i[I_i[i] + j];
260 
261  J[I[i+nrows] + j] = J_i[I_i[i] + j];
262  D[I[i+nrows] + j] = factor*D_i[I_i[i] + j];
263  }
264  }
265  }
266 
267  return new SparseMatrix(I, J, D, this->Height(), this->Width());
268 }
269 
270 
271 #ifdef MFEM_USE_SUITESPARSE
272 
274 {
275  mat = NULL;
276  Numeric = NULL;
277  AI = AJ = NULL;
278  if (!use_long_ints)
279  {
280  umfpack_zi_defaults(Control);
281  }
282  else
283  {
284  umfpack_zl_defaults(Control);
285  }
286 }
287 
289 {
290  int *Ap, *Ai;
291  void *Symbolic;
292  double *Ax;
293  double *Az;
294 
295  if (Numeric)
296  {
297  if (!use_long_ints)
298  {
299  umfpack_zi_free_numeric(&Numeric);
300  }
301  else
302  {
303  umfpack_zl_free_numeric(&Numeric);
304  }
305  }
306 
307  mat = const_cast<ComplexSparseMatrix *>
308  (dynamic_cast<const ComplexSparseMatrix *>(&op));
309  MFEM_VERIFY(mat, "not a ComplexSparseMatrix");
310 
311  MFEM_VERIFY(mat->real().NumNonZeroElems() == mat->imag().NumNonZeroElems(),
312  "Real and imag Sparsity pattern mismatch: Try setting Assemble (skip_zeros = 0)");
313 
314  // UMFPack requires that the column-indices in mat corresponding to each
315  // row be sorted.
316  // Generally, this will modify the ordering of the entries of mat.
317 
320 
321  height = mat->real().Height();
322  width = mat->real().Width();
323  MFEM_VERIFY(width == height, "not a square matrix");
324 
325  Ap = mat->real().GetI(); // assuming real and imag have the same sparsity
326  Ai = mat->real().GetJ();
327  Ax = mat->real().GetData();
328  Az = mat->imag().GetData();
329 
330  if (!use_long_ints)
331  {
332  int status = umfpack_zi_symbolic(width,width,Ap,Ai,Ax,Az,&Symbolic,
333  Control,Info);
334  if (status < 0)
335  {
336  umfpack_zi_report_info(Control, Info);
337  umfpack_zi_report_status(Control, status);
338  mfem_error("ComplexUMFPackSolver::SetOperator :"
339  " umfpack_zi_symbolic() failed!");
340  }
341 
342  status = umfpack_zi_numeric(Ap, Ai, Ax, Az, Symbolic, &Numeric,
343  Control, Info);
344  if (status < 0)
345  {
346  umfpack_zi_report_info(Control, Info);
347  umfpack_zi_report_status(Control, status);
348  mfem_error("ComplexUMFPackSolver::SetOperator :"
349  " umfpack_zi_numeric() failed!");
350  }
351  umfpack_zi_free_symbolic(&Symbolic);
352  }
353  else
354  {
355  SuiteSparse_long status;
356 
357  delete [] AJ;
358  delete [] AI;
359  AI = new SuiteSparse_long[width + 1];
360  AJ = new SuiteSparse_long[Ap[width]];
361  for (int i = 0; i <= width; i++)
362  {
363  AI[i] = (SuiteSparse_long)(Ap[i]);
364  }
365  for (int i = 0; i < Ap[width]; i++)
366  {
367  AJ[i] = (SuiteSparse_long)(Ai[i]);
368  }
369 
370  status = umfpack_zl_symbolic(width, width, AI, AJ, Ax, Az, &Symbolic,
371  Control, Info);
372  if (status < 0)
373  {
374  umfpack_zl_report_info(Control, Info);
375  umfpack_zl_report_status(Control, status);
376  mfem_error("ComplexUMFPackSolver::SetOperator :"
377  " umfpack_zl_symbolic() failed!");
378  }
379 
380  status = umfpack_zl_numeric(AI, AJ, Ax, Az, Symbolic, &Numeric,
381  Control, Info);
382  if (status < 0)
383  {
384  umfpack_zl_report_info(Control, Info);
385  umfpack_zl_report_status(Control, status);
386  mfem_error("ComplexUMFPackSolver::SetOperator :"
387  " umfpack_zl_numeric() failed!");
388  }
389  umfpack_zl_free_symbolic(&Symbolic);
390  }
391 }
392 
394 {
395  if (mat == NULL)
396  mfem_error("ComplexUMFPackSolver::Mult : matrix is not set!"
397  " Call SetOperator first!");
398  int n = b.Size()/2;
399  double * datax = x.GetData();
400  double * datab = b.GetData();
401 
402  // For the Block Symmetric case data the imaginary part
403  // has to be scaled by -1
405  Vector bimag;
406  if (conv == ComplexOperator::Convention::BLOCK_SYMMETRIC)
407  {
408  bimag.SetDataAndSize(&datab[n],n);
409  bimag *=-1.0;
410  }
411 
412  // Solve the transpose, since UMFPack expects CCS instead of CRS format
413  if (!use_long_ints)
414  {
415  int status =
416  umfpack_zi_solve(UMFPACK_Aat, mat->real().GetI(), mat->real().GetJ(),
417  mat->real().GetData(), mat->imag().GetData(),
418  datax, &datax[n], datab, &datab[n], Numeric, Control, Info);
419  umfpack_zi_report_info(Control, Info);
420  if (status < 0)
421  {
422  umfpack_zi_report_status(Control, status);
423  mfem_error("ComplexUMFPackSolver::Mult : umfpack_zi_solve() failed!");
424  }
425  }
426  else
427  {
428  SuiteSparse_long status =
429  umfpack_zl_solve(UMFPACK_Aat,AI,AJ,mat->real().GetData(),
430  mat->imag().GetData(),
431  datax,&datax[n],datab,&datab[n],Numeric,Control,Info);
432 
433  umfpack_zl_report_info(Control, Info);
434  if (status < 0)
435  {
436  umfpack_zl_report_status(Control, status);
437  mfem_error("ComplexUMFPackSolver::Mult : umfpack_zl_solve() failed!");
438  }
439  }
440  if (conv == ComplexOperator::Convention::BLOCK_SYMMETRIC)
441  {
442  bimag *=-1.0;
443  }
444 }
445 
447 {
448  if (mat == NULL)
449  mfem_error("ComplexUMFPackSolver::Mult : matrix is not set!"
450  " Call SetOperator first!");
451  int n = b.Size()/2;
452  double * datax = x.GetData();
453  double * datab = b.GetData();
454 
456  Vector bimag;
457  bimag.SetDataAndSize(&datab[n],n);
458 
459  // Solve the Adjoint A^H x = b by solving
460  // the conjugate problem A^T \bar{x} = \bar{b}
461  if ((!transa && conv == ComplexOperator::HERMITIAN) ||
463  {
464  bimag *=-1.0;
465  }
466 
467  if (!use_long_ints)
468  {
469  int status =
470  umfpack_zi_solve(UMFPACK_A, mat->real().GetI(), mat->real().GetJ(),
471  mat->real().GetData(), mat->imag().GetData(),
472  datax, &datax[n], datab, &datab[n], Numeric, Control, Info);
473  umfpack_zi_report_info(Control, Info);
474  if (status < 0)
475  {
476  umfpack_zi_report_status(Control, status);
477  mfem_error("ComplexUMFPackSolver::Mult : umfpack_zi_solve() failed!");
478  }
479  }
480  else
481  {
482  SuiteSparse_long status =
483  umfpack_zl_solve(UMFPACK_A,AI,AJ,mat->real().GetData(),
484  mat->imag().GetData(),
485  datax,&datax[n],datab,&datab[n],Numeric,Control,Info);
486 
487  umfpack_zl_report_info(Control, Info);
488  if (status < 0)
489  {
490  umfpack_zl_report_status(Control, status);
491  mfem_error("ComplexUMFPackSolver::Mult : umfpack_zl_solve() failed!");
492  }
493  }
494  if (!transa)
495  {
496  Vector ximag;
497  ximag.SetDataAndSize(&datax[n],n);
498  ximag *=-1.0;
499  }
500  if ((!transa && conv == ComplexOperator::HERMITIAN) ||
502  {
503  bimag *=-1.0;
504  }
505 }
506 
508 {
509  delete [] AJ;
510  delete [] AI;
511  if (Numeric)
512  {
513  if (!use_long_ints)
514  {
515  umfpack_zi_free_numeric(&Numeric);
516  }
517  else
518  {
519  umfpack_zl_free_numeric(&Numeric);
520  }
521  }
522 }
523 
524 #endif
525 
526 #ifdef MFEM_USE_MPI
527 
529  HypreParMatrix * A_Imag,
530  bool ownReal, bool ownImag,
531  Convention convention)
532  : ComplexOperator(A_Real, A_Imag, ownReal, ownImag, convention)
533 {
534  comm_ = (A_Real) ? A_Real->GetComm() :
535  ((A_Imag) ? A_Imag->GetComm() : MPI_COMM_WORLD);
536 
537  MPI_Comm_rank(comm_, &myid_);
538  MPI_Comm_size(comm_, &nranks_);
539 }
540 
542 {
543  MFEM_ASSERT(Op_Real_, "ComplexHypreParMatrix has no real part!");
544  return dynamic_cast<HypreParMatrix &>(*Op_Real_);
545 }
546 
548 {
549  MFEM_ASSERT(Op_Imag_, "ComplexHypreParMatrix has no imaginary part!");
550  return dynamic_cast<HypreParMatrix &>(*Op_Imag_);
551 }
552 
554 {
555  MFEM_ASSERT(Op_Real_, "ComplexHypreParMatrix has no real part!");
556  return dynamic_cast<const HypreParMatrix &>(*Op_Real_);
557 }
558 
560 {
561  MFEM_ASSERT(Op_Imag_, "ComplexHypreParMatrix has no imaginary part!");
562  return dynamic_cast<const HypreParMatrix &>(*Op_Imag_);
563 }
564 
566 {
567  HypreParMatrix * A_r = dynamic_cast<HypreParMatrix*>(Op_Real_);
568  HypreParMatrix * A_i = dynamic_cast<HypreParMatrix*>(Op_Imag_);
569 
570  if ( A_r == NULL && A_i == NULL ) { return NULL; }
571 
572  HYPRE_Int global_num_rows_r = (A_r) ? A_r->GetGlobalNumRows() : 0;
573  HYPRE_Int global_num_rows_i = (A_i) ? A_i->GetGlobalNumRows() : 0;
574  HYPRE_Int global_num_rows = std::max(global_num_rows_r, global_num_rows_i);
575 
576  HYPRE_Int global_num_cols_r = (A_r) ? A_r->GetGlobalNumCols() : 0;
577  HYPRE_Int global_num_cols_i = (A_i) ? A_i->GetGlobalNumCols() : 0;
578  HYPRE_Int global_num_cols = std::max(global_num_cols_r, global_num_cols_i);
579 
580  int row_starts_size = (HYPRE_AssumedPartitionCheck()) ? 2 : nranks_ + 1;
581  HYPRE_Int * row_starts = mfem_hypre_CTAlloc(HYPRE_Int, row_starts_size);
582  HYPRE_Int * col_starts = mfem_hypre_CTAlloc(HYPRE_Int, row_starts_size);
583 
584  const HYPRE_Int * row_starts_z = (A_r) ? A_r->RowPart() :
585  ((A_i) ? A_i->RowPart() : NULL);
586  const HYPRE_Int * col_starts_z = (A_r) ? A_r->ColPart() :
587  ((A_i) ? A_i->ColPart() : NULL);
588 
589  for (int i = 0; i < row_starts_size; i++)
590  {
591  row_starts[i] = 2 * row_starts_z[i];
592  col_starts[i] = 2 * col_starts_z[i];
593  }
594 
595  SparseMatrix diag_r, diag_i, offd_r, offd_i;
596  HYPRE_Int * cmap_r, * cmap_i;
597 
598  int nrows_r = 0, nrows_i = 0, ncols_r = 0, ncols_i = 0;
599  int ncols_offd_r = 0, ncols_offd_i = 0;
600  if (A_r)
601  {
602  A_r->GetDiag(diag_r);
603  A_r->GetOffd(offd_r, cmap_r);
604  nrows_r = diag_r.Height();
605  ncols_r = diag_r.Width();
606  ncols_offd_r = offd_r.Width();
607  }
608  if (A_i)
609  {
610  A_i->GetDiag(diag_i);
611  A_i->GetOffd(offd_i, cmap_i);
612  nrows_i = diag_i.Height();
613  ncols_i = diag_i.Width();
614  ncols_offd_i = offd_i.Width();
615  }
616  int nrows = std::max(nrows_r, nrows_i);
617  int ncols = std::max(ncols_r, ncols_i);
618 
619  // Determine the unique set of off-diagonal columns global indices
620  std::set<int> cset;
621  for (int i=0; i<ncols_offd_r; i++)
622  {
623  cset.insert(cmap_r[i]);
624  }
625  for (int i=0; i<ncols_offd_i; i++)
626  {
627  cset.insert(cmap_i[i]);
628  }
629  int num_cols_offd = (int)cset.size();
630 
631  // Extract pointers to the various CSR arrays of the diagonal blocks
632  const int * diag_r_I = (A_r) ? diag_r.GetI() : NULL;
633  const int * diag_i_I = (A_i) ? diag_i.GetI() : NULL;
634 
635  const int * diag_r_J = (A_r) ? diag_r.GetJ() : NULL;
636  const int * diag_i_J = (A_i) ? diag_i.GetJ() : NULL;
637 
638  const double * diag_r_D = (A_r) ? diag_r.GetData() : NULL;
639  const double * diag_i_D = (A_i) ? diag_i.GetData() : NULL;
640 
641  int diag_r_nnz = (diag_r_I) ? diag_r_I[nrows] : 0;
642  int diag_i_nnz = (diag_i_I) ? diag_i_I[nrows] : 0;
643  int diag_nnz = 2 * (diag_r_nnz + diag_i_nnz);
644 
645  // Extract pointers to the various CSR arrays of the off-diagonal blocks
646  const int * offd_r_I = (A_r) ? offd_r.GetI() : NULL;
647  const int * offd_i_I = (A_i) ? offd_i.GetI() : NULL;
648 
649  const int * offd_r_J = (A_r) ? offd_r.GetJ() : NULL;
650  const int * offd_i_J = (A_i) ? offd_i.GetJ() : NULL;
651 
652  const double * offd_r_D = (A_r) ? offd_r.GetData() : NULL;
653  const double * offd_i_D = (A_i) ? offd_i.GetData() : NULL;
654 
655  int offd_r_nnz = (offd_r_I) ? offd_r_I[nrows] : 0;
656  int offd_i_nnz = (offd_i_I) ? offd_i_I[nrows] : 0;
657  int offd_nnz = 2 * (offd_r_nnz + offd_i_nnz);
658 
659  // Allocate CSR arrays for the combined matrix
660  HYPRE_Int * diag_I = mfem_hypre_CTAlloc(HYPRE_Int, 2 * nrows + 1);
661  HYPRE_Int * diag_J = mfem_hypre_CTAlloc(HYPRE_Int, diag_nnz);
662  double * diag_D = mfem_hypre_CTAlloc(double, diag_nnz);
663 
664  HYPRE_Int * offd_I = mfem_hypre_CTAlloc(HYPRE_Int, 2 * nrows + 1);
665  HYPRE_Int * offd_J = mfem_hypre_CTAlloc(HYPRE_Int, offd_nnz);
666  double * offd_D = mfem_hypre_CTAlloc(double, offd_nnz);
667  HYPRE_Int * cmap = mfem_hypre_CTAlloc(HYPRE_Int, 2 * num_cols_offd);
668 
669  // Fill the CSR arrays for the diagonal portion of the matrix
670  const double factor = (convention_ == HERMITIAN) ? 1.0 : -1.0;
671 
672  diag_I[0] = 0;
673  diag_I[nrows] = diag_r_nnz + diag_i_nnz;
674  for (int i=0; i<nrows; i++)
675  {
676  diag_I[i + 1] = ((diag_r_I)?diag_r_I[i+1]:0) +
677  ((diag_i_I)?diag_i_I[i+1]:0);
678  diag_I[i + nrows + 1] = diag_I[i+1] + diag_r_nnz + diag_i_nnz;
679 
680  if (diag_r_I)
681  {
682  for (int j=0; j<diag_r_I[i+1] - diag_r_I[i]; j++)
683  {
684  diag_J[diag_I[i] + j] = diag_r_J[diag_r_I[i] + j];
685  diag_D[diag_I[i] + j] = diag_r_D[diag_r_I[i] + j];
686 
687  diag_J[diag_I[i+nrows] + j] =
688  diag_r_J[diag_r_I[i] + j] + ncols;
689  diag_D[diag_I[i+nrows] + j] =
690  factor * diag_r_D[diag_r_I[i] + j];
691  }
692  }
693  if (diag_i_I)
694  {
695  const int off_r = (diag_r_I)?(diag_r_I[i+1] - diag_r_I[i]):0;
696  for (int j=0; j<diag_i_I[i+1] - diag_i_I[i]; j++)
697  {
698  diag_J[diag_I[i] + off_r + j] = diag_i_J[diag_i_I[i] + j] + ncols;
699  diag_D[diag_I[i] + off_r + j] = -diag_i_D[diag_i_I[i] + j];
700 
701  diag_J[diag_I[i+nrows] + off_r + j] = diag_i_J[diag_i_I[i] + j];
702  diag_D[diag_I[i+nrows] + off_r + j] =
703  factor * diag_i_D[diag_i_I[i] + j];
704  }
705  }
706  }
707 
708  // Determine the mappings describing the layout of off-diagonal columns
709  int num_recv_procs = 0;
710  HYPRE_Int * offd_col_start_stop = NULL;
711  this->getColStartStop(A_r, A_i, num_recv_procs, offd_col_start_stop);
712 
713  std::set<int>::iterator sit;
714  std::map<int,int> cmapa, cmapb, cinvmap;
715  for (sit=cset.begin(); sit!=cset.end(); sit++)
716  {
717  int col_orig = *sit;
718  int col_2x2 = -1;
719  int col_size = 0;
720  for (int i=0; i<num_recv_procs; i++)
721  {
722  if (offd_col_start_stop[2*i] <= col_orig &&
723  col_orig < offd_col_start_stop[2*i+1])
724  {
725  col_2x2 = offd_col_start_stop[2*i] + col_orig;
726  col_size = offd_col_start_stop[2*i+1] - offd_col_start_stop[2*i];
727  break;
728  }
729  }
730  cmapa[*sit] = col_2x2;
731  cmapb[*sit] = col_2x2 + col_size;
732  cinvmap[col_2x2] = -1;
733  cinvmap[col_2x2 + col_size] = -1;
734  }
735  delete [] offd_col_start_stop;
736 
737  std::map<int, int>::iterator mit;
738  int i = 0;
739  for (mit=cinvmap.begin(); mit!=cinvmap.end(); mit++, i++)
740  {
741  mit->second = i;
742  cmap[i] = mit->first;
743  }
744 
745  // Fill the CSR arrays for the off-diagonal portion of the matrix
746  offd_I[0] = 0;
747  offd_I[nrows] = offd_r_nnz + offd_i_nnz;
748  for (int i=0; i<nrows; i++)
749  {
750  offd_I[i + 1] = ((offd_r_I)?offd_r_I[i+1]:0) +
751  ((offd_i_I)?offd_i_I[i+1]:0);
752  offd_I[i + nrows + 1] = offd_I[i+1] + offd_r_nnz + offd_i_nnz;
753 
754  if (offd_r_I)
755  {
756  const int off_i = (offd_i_I)?(offd_i_I[i+1] - offd_i_I[i]):0;
757  for (int j=0; j<offd_r_I[i+1] - offd_r_I[i]; j++)
758  {
759  offd_J[offd_I[i] + j] =
760  cinvmap[cmapa[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
761  offd_D[offd_I[i] + j] = offd_r_D[offd_r_I[i] + j];
762 
763  offd_J[offd_I[i+nrows] + off_i + j] =
764  cinvmap[cmapb[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
765  offd_D[offd_I[i+nrows] + off_i + j] =
766  factor * offd_r_D[offd_r_I[i] + j];
767  }
768  }
769  if (offd_i_I)
770  {
771  const int off_r = (offd_r_I)?(offd_r_I[i+1] - offd_r_I[i]):0;
772  for (int j=0; j<offd_i_I[i+1] - offd_i_I[i]; j++)
773  {
774  offd_J[offd_I[i] + off_r + j] =
775  cinvmap[cmapb[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
776  offd_D[offd_I[i] + off_r + j] = -offd_i_D[offd_i_I[i] + j];
777 
778  offd_J[offd_I[i+nrows] + j] =
779  cinvmap[cmapa[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
780  offd_D[offd_I[i+nrows] + j] = factor * offd_i_D[offd_i_I[i] + j];
781  }
782  }
783  }
784 
785  // Construct the combined matrix
786  HypreParMatrix * A = new HypreParMatrix(comm_,
787  2 * global_num_rows,
788  2 * global_num_cols,
789  row_starts, col_starts,
790  diag_I, diag_J, diag_D,
791  offd_I, offd_J, offd_D,
792  2 * num_cols_offd, cmap);
793 
794  // Give the new matrix ownership of its internal arrays
795  A->SetOwnerFlags(-1,-1,-1);
796  hypre_CSRMatrixSetDataOwner(((hypre_ParCSRMatrix*)(*A))->diag,1);
797  hypre_CSRMatrixSetDataOwner(((hypre_ParCSRMatrix*)(*A))->offd,1);
798  hypre_ParCSRMatrixSetRowStartsOwner((hypre_ParCSRMatrix*)(*A),1);
799  hypre_ParCSRMatrixSetColStartsOwner((hypre_ParCSRMatrix*)(*A),1);
800 
801  return A;
802 }
803 
804 void
805 ComplexHypreParMatrix::getColStartStop(const HypreParMatrix * A_r,
806  const HypreParMatrix * A_i,
807  int & num_recv_procs,
808  HYPRE_Int *& offd_col_start_stop) const
809 {
810  hypre_ParCSRCommPkg * comm_pkg_r =
811  (A_r) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_r)) : NULL;
812  hypre_ParCSRCommPkg * comm_pkg_i =
813  (A_i) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_i)) : NULL;
814 
815  std::set<HYPRE_Int> send_procs, recv_procs;
816  if ( comm_pkg_r )
817  {
818  for (HYPRE_Int i=0; i<comm_pkg_r->num_sends; i++)
819  {
820  send_procs.insert(comm_pkg_r->send_procs[i]);
821  }
822  for (HYPRE_Int i=0; i<comm_pkg_r->num_recvs; i++)
823  {
824  recv_procs.insert(comm_pkg_r->recv_procs[i]);
825  }
826  }
827  if ( comm_pkg_i )
828  {
829  for (HYPRE_Int i=0; i<comm_pkg_i->num_sends; i++)
830  {
831  send_procs.insert(comm_pkg_i->send_procs[i]);
832  }
833  for (HYPRE_Int i=0; i<comm_pkg_i->num_recvs; i++)
834  {
835  recv_procs.insert(comm_pkg_i->recv_procs[i]);
836  }
837  }
838 
839  num_recv_procs = (int)recv_procs.size();
840 
841  HYPRE_Int loc_start_stop[2];
842  offd_col_start_stop = new HYPRE_Int[2 * num_recv_procs];
843 
844  const HYPRE_Int * row_part = (A_r) ? A_r->RowPart() :
845  ((A_i) ? A_i->RowPart() : NULL);
846 
847  int row_part_ind = (HYPRE_AssumedPartitionCheck()) ? 0 : myid_;
848  loc_start_stop[0] = row_part[row_part_ind];
849  loc_start_stop[1] = row_part[row_part_ind+1];
850 
851  MPI_Request * req = new MPI_Request[send_procs.size()+recv_procs.size()];
852  MPI_Status * stat = new MPI_Status[send_procs.size()+recv_procs.size()];
853  int send_count = 0;
854  int recv_count = 0;
855  int tag = 0;
856 
857  std::set<HYPRE_Int>::iterator sit;
858  for (sit=send_procs.begin(); sit!=send_procs.end(); sit++)
859  {
860  MPI_Isend(loc_start_stop, 2, HYPRE_MPI_INT,
861  *sit, tag, comm_, &req[send_count]);
862  send_count++;
863  }
864  for (sit=recv_procs.begin(); sit!=recv_procs.end(); sit++)
865  {
866  MPI_Irecv(&offd_col_start_stop[2*recv_count], 2, HYPRE_MPI_INT,
867  *sit, tag, comm_, &req[send_count+recv_count]);
868  recv_count++;
869  }
870 
871  MPI_Waitall(send_count+recv_count, req, stat);
872 
873  delete [] req;
874  delete [] stat;
875 }
876 
877 #endif // MFEM_USE_MPI
878 
879 }
double Control[UMFPACK_CONTROL]
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1359
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:330
virtual Operator & real()
Real or imaginary part accessor methods.
ComplexSparseMatrix * mat
ComplexHypreParMatrix(HypreParMatrix *A_Real, HypreParMatrix *A_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:370
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
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:71
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:178
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:89
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
int * GetI()
Return the array I.
Definition: sparsemat.hpp:173
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:169
ComplexOperator(Operator *Op_Real, Operator *Op_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
Constructs complex operator object.
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:409
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:92
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:183
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:419
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
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:342
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:153
double b
Definition: lissajous.cpp:42
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:887
virtual SparseMatrix & real()
Real or imaginary part accessor methods.
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:366
SparseMatrix * GetSystemMatrix() const
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
virtual Operator & imag()
Convention GetConvention() const
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:128
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a ComplexSparseMatrix.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:213
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
virtual void Mult(const Vector &b, Vector &x) const
This is solving the system A x = b.
virtual SparseMatrix & imag()
void Destroy()
Destroy a vector.
Definition: vector.hpp:530
Native convention for Hermitian operators.
virtual void MultTranspose(const Vector &b, Vector &x) const
This is solving the system: A^H x = b (when transa = false) This is equivalent to solving the transpo...
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:51
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:517
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:923
HypreParMatrix * GetSystemMatrix() const
virtual HypreParMatrix & imag()
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
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
void SyncAliasMemory(const Vector &v)
Update the alias memory location of the vector to match v.
Definition: vector.hpp:193