MFEM  v3.4
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, 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 "complex_operator.hpp"
13 
14 namespace mfem
15 {
16 
18  bool ownReal, bool ownImag,
19  Convention convention)
20  : Operator(2*Op_Real->Height(), 2*Op_Real->Width())
21  , Op_Real_(Op_Real)
22  , Op_Imag_(Op_Imag)
23  , ownReal_(ownReal)
24  , ownImag_(ownImag)
25  , convention_(convention)
26  , x_r_(NULL, Op_Real->Width())
27  , x_i_(NULL, Op_Real->Width())
28  , y_r_(NULL, Op_Real->Height())
29  , y_i_(NULL, Op_Real->Height())
30  , u_(NULL)
31  , v_(NULL)
32 {}
33 
35 {
36  if (ownReal_) { delete Op_Real_; }
37  if (ownImag_) { delete Op_Imag_; }
38  delete u_;
39  delete v_;
40 }
41 
42 void ComplexOperator::Mult(const Vector &x, Vector &y) const
43 {
44  double * x_data = x.GetData();
45  x_r_.SetData(x_data);
46  x_i_.SetData(&x_data[Op_Real_->Width()]);
47 
48  y_r_.SetData(&y[0]);
49  y_i_.SetData(&y[Op_Real_->Height()]);
50 
51  this->Mult(x_r_, x_i_, y_r_, y_i_);
52 }
53 
54 void ComplexOperator::Mult(const Vector &x_r, const Vector &x_i,
55  Vector &y_r, Vector &y_i) const
56 {
57  if (Op_Real_)
58  {
59  Op_Real_->Mult(x_r, y_r);
60  Op_Real_->Mult(x_i, y_i);
61  }
62  else
63  {
64  y_r = 0.0;
65  y_i = 0.0;
66  }
67  if (Op_Imag_)
68  {
69  if (!v_) { v_ = new Vector(Op_Imag_->Height()); }
70  Op_Imag_->Mult(x_i, *v_);
71  y_r_ -= *v_;
72  Op_Imag_->Mult(x_r, *v_);
73  y_i_ += *v_;
74  }
75 
77  {
78  y_i_ *= -1.0;
79  }
80 }
81 
83 {
84  double * x_data = x.GetData();
85  y_r_.SetData(x_data);
86  y_i_.SetData(&x_data[Op_Real_->Height()]);
87 
88  x_r_.SetData(&y[0]);
89  x_i_.SetData(&y[Op_Real_->Width()]);
90 
91  this->MultTranspose(y_r_, y_i_, x_r_, x_i_);
92 }
93 
94 void ComplexOperator::MultTranspose(const Vector &x_r, const Vector &x_i,
95  Vector &y_r, Vector &y_i) const
96 {
97  if (Op_Real_)
98  {
99  Op_Real_->MultTranspose(x_r, y_r);
100  Op_Real_->MultTranspose(x_i, y_i);
101 
103  {
104  y_i *= -1.0;
105  }
106  }
107  else
108  {
109  y_r = 0.0;
110  y_i = 0.0;
111  }
112  if (Op_Imag_)
113  {
114  if (!u_) { u_ = new Vector(Op_Imag_->Width()); }
115  Op_Imag_->MultTranspose(x_i, *u_);
116  y_r_.Add(convention_ == BLOCK_SYMMETRIC ? -1.0 : 1.0, *u_);
117  Op_Imag_->MultTranspose(x_r, *u_);
118  y_i_ -= *u_;
119  }
120 }
121 
122 
124 {
125  SparseMatrix * A_r = dynamic_cast<SparseMatrix*>(Op_Real_);
126  SparseMatrix * A_i = dynamic_cast<SparseMatrix*>(Op_Imag_);
127 
128  const int nrows_r = (A_r)?A_r->Height():0;
129  const int nrows_i = (A_i)?A_i->Height():0;
130  const int nrows = std::max(nrows_r, nrows_i);
131 
132  const int *I_r = (A_r)?A_r->GetI():NULL;
133  const int *I_i = (A_i)?A_i->GetI():NULL;
134 
135  const int *J_r = (A_r)?A_r->GetJ():NULL;
136  const int *J_i = (A_i)?A_i->GetJ():NULL;
137 
138  const double *D_r = (A_r)?A_r->GetData():NULL;
139  const double *D_i = (A_i)?A_i->GetData():NULL;
140 
141  const int nnz_r = (I_r)?I_r[nrows]:0;
142  const int nnz_i = (I_i)?I_i[nrows]:0;
143  const int nnz = 2 * (nnz_r + nnz_i);
144 
145  int *I = new int[this->Height()+1];
146  int *J = new int[nnz];
147  double *D = new double[nnz];
148 
149  const double factor = (convention_ == HERMITIAN) ? 1.0 : -1.0;
150 
151  I[0] = 0;
152  I[nrows] = nnz_r + nnz_i;
153  for (int i=0; i<nrows; i++)
154  {
155  I[i + 1] = ((I_r)?I_r[i+1]:0) + ((I_i)?I_i[i+1]:0);
156  I[i + nrows + 1] = I[i+1] + nnz_r + nnz_i;
157 
158  if (I_r)
159  {
160  const int off_i = (I_i)?(I_i[i+1] - I_i[i]):0;
161  for (int j=0; j<I_r[i+1] - I_r[i]; j++)
162  {
163  J[I[i] + j] = J_r[I_r[i] + j];
164  D[I[i] + j] = D_r[I_r[i] + j];
165 
166  J[I[i+nrows] + off_i + j] = J_r[I_r[i] + j] + nrows;
167  D[I[i+nrows] + off_i + j] = factor*D_r[I_r[i] + j];
168  }
169  }
170  if (I_i)
171  {
172  const int off_r = (I_r)?(I_r[i+1] - I_r[i]):0;
173  for (int j=0; j<I_i[i+1] - I_i[i]; j++)
174  {
175  J[I[i] + off_r + j] = J_i[I_i[i] + j] + nrows;
176  D[I[i] + off_r + j] = -D_i[I_i[i] + j];
177 
178  J[I[i+nrows] + j] = J_i[I_i[i] + j];
179  D[I[i+nrows] + j] = factor*D_i[I_i[i] + j];
180  }
181  }
182  }
183 
184  return new SparseMatrix(I, J, D, this->Height(), this->Width());
185 }
186 
187 }
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:42
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:129
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:52
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
void SetData(double *d)
Definition: vector.hpp:94
double * GetData() const
Return element data, i.e. array A.
Definition: sparsemat.hpp:138
int * GetI() const
Return the array I.
Definition: sparsemat.hpp:134
SparseMatrix * GetSystemMatrix() const
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:201
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
Abstract operator.
Definition: operator.hpp:21
int * GetJ() const
Return the array J.
Definition: sparsemat.hpp:136