MAST
lapack_zggev_interface.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 // MAST includes
23 
24 
25 void
27  const ComplexMatrixX &B,
28  bool computeEigenvectors) {
29 
30  libmesh_assert(A.cols() == A.rows() &&
31  B.cols() == A.rows() &&
32  B.cols() == B.rows());
33 
34  _A = A;
35  _B = B;
36 
38  Amat = _A,
39  Bmat = _B;
40 
41  int n = (int)A.cols();
42 
43  char L='N',R='N';
44 
45  if (computeEigenvectors)
46  {
47  L = 'V'; R = 'V';
48  VL.setZero(n, n);
49  VR.setZero(n, n);
50  }
51 
52  int
53  lwork=16*n,
54  l_rwork=8*n;
55  info_val=-1;
56 
57  alpha.setZero(n);
58  beta.setZero(n);
60  work = ComplexVectorX::Zero(lwork);
62  rwork = RealVectorX::Zero(l_rwork);
63 
64  Complex
65  *a_vals = Amat.data(),
66  *b_vals = Bmat.data(),
67  *alpha_v = alpha.data(),
68  *beta_v = beta.data(),
69  *VL_v = VL.data(),
70  *VR_v = VR.data(),
71  *work_v = work.data();
72 
73  Real
74  *rwork_v = rwork.data();
75 
76 
77  zggev_(&L, &R, &n,
78  &(a_vals[0]), &n,
79  &(b_vals[0]), &n,
80  &(alpha_v[0]), &(beta_v[0]),
81  &(VL_v[0]), &n, &(VR_v[0]), &n,
82  &(work_v[0]), &lwork,
83  &(rwork_v[0]),
84  &info_val);
85 
86  if (info_val != 0)
87  libMesh::out
88  << "Warning!! ZGGEV returned with nonzero info = "
89  << info_val << std::endl;
90 }
91 
92 
virtual void compute(const ComplexMatrixX &A, const ComplexMatrixX &B, bool computeEigenvectors=true)
computes the eigensolution for .
Matrix< Complex, Dynamic, 1 > ComplexVectorX
const ComplexMatrixX & B() const
libMesh::Real Real
libMesh::Complex Complex
Matrix< Complex, Dynamic, Dynamic > ComplexMatrixX
Matrix< Real, Dynamic, 1 > RealVectorX
const ComplexMatrixX & A() const
int zggev_(char *jobvl, char *jobvr, int *n, std::complex< double > *a, int *lda, std::complex< double > *b, int *ldb, std::complex< double > *alpha, std::complex< double > *beta, std::complex< double > *vl, int *ldvl, std::complex< double > *vr, int *ldvr, std::complex< double > *work, int *lwork, double *rwork, int *info)