1 #ifndef SELDON_FILE_FEAST_CXX
2 #define SELDON_FILE_FEAST_CXX
8 #ifdef SELDON_WITH_VIRTUAL
10 void FindEigenvaluesFeast(EigenProblem_Base<T>& var,
11 Vector<T>& eigen_values,
12 Vector<T>& lambda_imag,
13 Matrix<T, General, ColMajor>& eigen_vectors)
15 template<
class EigenProblem,
class T,
class Allocator1,
16 class Allocator2,
class Allocator3>
17 void FindEigenvaluesFeast(EigenProblem& var,
18 Vector<T, VectFull, Allocator1>& eigen_values,
19 Vector<T, VectFull, Allocator2>& lambda_imag,
20 Matrix<T, General, ColMajor, Allocator3>& eigen_vectors)
26 feastinit_(fpm.GetData());
29 if (var.GetPrintLevel() > 0)
32 double tol = var.GetStoppingCriterion();
35 for (
int i = 64; i < fpm.GetM(); i++)
38 fpm(3) = var.GetNbMaximumIterations();
40 FeastParam& param = var.GetFeastParameters();
41 if (param.EstimateNumberEigenval())
46 var.SetComputationalMode(var.REGULAR_MODE);
48 #ifdef SELDON_WITH_MPI
50 MPI_Comm_dup(var.GetCommunicator(), &solver_comm);
52 fpm(8) = MPI_Comm_c2f(var.GetRciCommunicator());
53 fpm(48) = MPI_Comm_c2f(solver_comm);
54 fpm(58) = MPI_Comm_c2f(var.GetGlobalCommunicator());
60 typedef typename ClassComplexType<T>::Tcplx Tcplx;
61 typedef typename ClassComplexType<T>::Treal Treal;
63 int m0 = var.GetNbAskedEigenvalues()+1;
65 Treal emin = param.GetLowerBoundInterval();
66 Treal emax = param.GetUpperBoundInterval();
67 Tcplx Emid = param.GetCircleCenterSpectrum();
68 Treal r = param.GetCircleRadiusSpectrum();
69 fpm(17) = param.GetRatioEllipseSpectrum();
70 fpm(18) = param.GetAngleEllipseSpectrum();
76 bool herm = var.IsHermitianProblem();
77 bool sym = var.IsSymmetricProblem();
80 if (var.UseCholeskyFactoForMass())
82 cout <<
"Cholesky case not implemented in Feast interface" << endl;
94 if (param.GetTypeIntegration() >= 0)
95 fpm(15) = param.GetTypeIntegration();
97 if (param.GetNumOfQuadraturePoints() > 0)
98 fpm(1) = param.GetNumOfQuadraturePoints();
104 if (param.GetTypeIntegration() >= 0)
105 fpm(15) = param.GetTypeIntegration();
107 if (param.GetNumOfQuadraturePoints() > 0)
108 fpm(7) = param.GetNumOfQuadraturePoints();
114 Matrix<T, General, ColMajor> work;
115 Matrix<Tcplx, General, ColMajor> zwork;
116 Matrix<Tcplx, General, ColMajor> workc;
117 Vector<Tcplx> xc(n), yc(n);
119 Vector<Tcplx> zaq, zbq;
120 Vector<T> x(n), y(n);
122 workc.Reallocate(n, m0);
125 work.Reallocate(n, m0);
126 aq.Reallocate(m0*m0);
127 bq.Reallocate(m0*m0);
131 zwork.Reallocate(n, m0);
132 zaq.Reallocate(m0*m0);
133 zbq.Reallocate(m0*m0);
137 zwork.Reallocate(n, m0);
138 zaq.Reallocate(m0*m0);
139 zbq.Reallocate(m0*m0);
142 work.Zero(); zwork.Zero();
143 workc.Zero(); xc.Zero(); yc.Zero();
145 aq.Zero(); bq.Zero(); zaq.Zero(); zbq.Zero();
148 Treal epsout(0); Tcplx ze;
149 Vector<Treal> res(m0);
152 Matrix<Tcplx, General, ColMajor> eigen_vec_cplx;
153 Vector<Treal> lambda; Vector<Tcplx> lambda_cplx;
157 eigen_vectors.Reallocate(n, m0);
158 lambda.Reallocate(m0);
162 eigen_vec_cplx.Reallocate(n, m0);
163 lambda_cplx.Reallocate(m0);
167 eigen_vec_cplx.Reallocate(n, 2*m0);
168 lambda_cplx.Reallocate(m0);
169 res.Reallocate(2*m0);
173 eigen_vectors.Zero(); lambda.Zero();
174 lambda_cplx.Zero(); eigen_vec_cplx.Zero();
176 if (var.DiagonalMass())
180 var.ComputeDiagonalMass();
183 var.FactorizeDiagonalMass();
185 else if (var.UseCholeskyFactoForMass())
190 var.ComputeMassForCholesky();
193 var.FactorizeCholeskyMass();
196 var.ComputeMassMatrix();
199 var.ComputeStiffnessMatrix();
202 int ijob = -1, info = 0;
206 CallFeast(ijob, n, herm, sym, ze, work, workc, aq, bq, zwork, zaq, zbq,
207 fpm, epsout, loop, emin, emax, Emid, r, m0, lambda, eigen_vectors,
208 lambda_cplx, eigen_vec_cplx, m, res, info);
213 var.ComputeAndFactorizeStiffnessMatrix(ze, -one);
218 for (
int k = 0; k < fpm(22); k++)
220 for (
int i = 0; i < n; i++)
223 if (var.DiagonalMass())
224 var.MltSqrtDiagonalMass(xc);
225 else if (var.UseCholeskyFactoForMass())
226 var.MltCholeskyMass(SeldonNoTrans, xc);
228 var.ComputeSolution(xc, yc);
229 var.IncrementLinearSolves();
231 if (var.DiagonalMass())
232 var.MltSqrtDiagonalMass(yc);
233 else if (var.UseCholeskyFactoForMass())
235 var.MltCholeskyMass(SeldonTrans, yc);
236 var.IncrementLinearSolves();
239 for (
int i = 0; i < n; i++)
251 for (
int k = 0; k < fpm(22); k++)
253 for (
int i = 0; i < n; i++)
256 if (var.DiagonalMass())
257 var.MltSqrtDiagonalMass(xc);
258 else if (var.UseCholeskyFactoForMass())
259 var.MltCholeskyMass(SeldonNoTrans, xc);
262 var.ComputeSolution(SeldonTrans, xc, yc);
263 var.IncrementLinearSolves();
266 if (var.DiagonalMass())
267 var.MltSqrtDiagonalMass(yc);
268 else if (var.UseCholeskyFactoForMass())
270 var.MltCholeskyMass(SeldonTrans, yc);
271 var.IncrementLinearSolves();
274 for (
int i = 0; i < n; i++)
281 int i = fpm(23), j = fpm(23) + fpm(24)-1;
283 for (
int k = i; k <= j; k++)
285 GetCol(eigen_vectors, k-1, x);
287 if (var.DiagonalMass())
288 var.MltInvSqrtDiagonalMass(x);
289 else if (var.UseCholeskyFactoForMass())
290 var.SolveCholeskyMass(SeldonTrans, x);
292 var.MltStiffness(x, y);
293 var.IncrementProdMatVect();
295 if (var.DiagonalMass())
296 var.MltInvSqrtDiagonalMass(y);
297 else if (var.UseCholeskyFactoForMass())
299 var.SolveCholeskyMass(SeldonNoTrans, y);
300 var.IncrementLinearSolves();
303 for (
int p = 0; p < n; p++)
307 for (
int k = i; k <= j; k++)
309 GetCol(eigen_vec_cplx, k-1, xc);
311 if (var.DiagonalMass())
312 var.MltInvSqrtDiagonalMass(xc);
314 var.MltStiffness(xc, yc);
315 var.IncrementProdMatVect();
317 if (var.DiagonalMass())
318 var.MltInvSqrtDiagonalMass(yc);
320 SetCol(yc, k-1, zwork);
326 int i = fpm(33), j = fpm(33) + fpm(34)-1;
329 cout <<
"impossible" << endl;
333 for (
int k = i; k <= j; k++)
335 GetCol(eigen_vec_cplx, k-1, xc);
337 if (var.DiagonalMass())
338 var.MltInvSqrtDiagonalMass(xc);
340 var.MltStiffness(SeldonConjTrans, xc, yc);
341 var.IncrementProdMatVect();
343 if (var.DiagonalMass())
344 var.MltInvSqrtDiagonalMass(yc);
346 SetCol(yc, k-1, zwork);
352 int i = fpm(23), j = fpm(23) + fpm(24)-1;
354 for (
int k = i; k <= j; k++)
356 for (
int p = 0; p < n; p++)
357 x(p) = eigen_vectors(p, k-1);
359 if (var.DiagonalMass() || var.UseCholeskyFactoForMass())
364 var.IncrementProdMatVect();
367 for (
int p = 0; p < n; p++)
371 for (
int k = i; k <= j; k++)
373 GetCol(eigen_vec_cplx, k-1, xc);
375 if (var.DiagonalMass())
380 var.IncrementProdMatVect();
383 SetCol(yc, k-1, zwork);
389 int i = fpm(33), j = fpm(33) + fpm(34)-1;
391 for (
int k = i; k <= j; k++)
393 for (
int p = 0; p < n; p++)
394 x(p) = eigen_vectors(p, k-1);
396 if (var.DiagonalMass() || var.UseCholeskyFactoForMass())
400 var.MltMass(SeldonConjTrans, x, y);
401 var.IncrementProdMatVect();
404 for (
int p = 0; p < n; p++)
408 for (
int k = i; k <= j; k++)
410 GetCol(eigen_vec_cplx, k-1, xc);
412 if (var.DiagonalMass())
416 var.MltMass(SeldonConjTrans, xc, yc);
417 var.IncrementProdMatVect();
420 SetCol(yc, k-1, zwork);
429 cout <<
"Feast returns error code 3 " << endl;
430 cout <<
"It usually means that there are more than "
431 << m0-1 <<
" eigenvalues in the ";
434 cout <<
"interval [" << emin <<
", " << emax <<
"] " << endl;
436 cout <<
"disk of center " << Emid <<
" and radius" << r << endl;
438 cout <<
"You can try to change the ";
444 cout <<
"or increase the number of asked eigenvalues" << endl;
450 cout <<
"No eigenvalues have been found in the ";
452 cout <<
"interval [" << emin <<
", " << emax <<
"] " << endl;
454 cout <<
"disk of center " << Emid <<
" and radius" << r << endl;
456 cout <<
"You can try to change the ";
458 cout <<
"interval. " << endl;
460 cout <<
"disk. " << endl;
466 cout <<
"Feast did not converge, loop = " << loop << endl;
467 cout <<
"Maximum number of iterations = " << var.GetNbMaximumIterations() << endl;
468 cout <<
"Residual = " << res << endl;
472 if (var.GetGlobalRankCommunicator() == 0)
473 cout <<
"Estimated number of eigenvalues = " << m << endl;
479 if (var.GetGlobalRankCommunicator() == 0)
480 cout <<
"Warning : subspace is not biorthogonal" << endl;
492 cout <<
"An error occurred during feast eigenvalue solving " << endl;
493 cout <<
"info = " << info << endl;
498 if ((var.GetPrintLevel() > 0) && (var.GetGlobalRankCommunicator() == 0))
499 cout <<
"Feast converged in " << loop <<
" cycles" << endl;
501 work.Clear(); workc.Clear(); zwork.Clear();
502 aq.Clear(); bq.Clear(); zaq.Clear(); zbq.Clear();
504 eigen_values.Reallocate(m);
505 lambda_imag.Reallocate(m);
509 for (
int i = 0; i < m; i++)
510 eigen_values(i) = lambda(i);
512 eigen_vectors.Resize(n, m);
515 var.FillComplexEigenvectors(m, Emid, tol, lambda_cplx, eigen_vec_cplx,
516 eigen_values, lambda_imag, eigen_vectors);
519 lambda_cplx.Clear(); eigen_vec_cplx.Clear();
521 T shiftr = var.GetShiftValue(), shifti = var.GetImagShiftValue();
528 void CallFeast(
int& ijob,
int n,
bool herm,
bool sym, complex<double>& ze,
529 Matrix<double, General, ColMajor>& work,
530 Matrix<complex<double>, General, ColMajor>& workc,
531 Vector<double>& aq, Vector<double>& bq,
532 Matrix<complex<double>, General, ColMajor>& zwork,
533 Vector<complex<double> >& zaq, Vector<complex<double> >& zbq,
534 IVect& fpm,
double& epsout,
int& loop,
535 double emin,
double emax, complex<double> Emid,
double r,
536 int m0, Vector<double>& lambda,
537 Matrix<double, General, ColMajor>& eigen_vectors,
538 Vector<complex<double> >& lambda_cplx,
539 Matrix<complex<double>, General, ColMajor>& eigen_vec_cplx,
540 int& m, Vector<double>& res,
int& info)
543 dfeast_srci_(&ijob, &n,
reinterpret_cast<void*
>(&ze), work.GetData(),
544 workc.GetDataVoid(), aq.GetData(), bq.GetData(),
545 fpm.GetData(), &epsout, &loop, &emin, &emax,
546 &m0, lambda.GetData(), eigen_vectors.GetData(),
547 &m, res.GetData(), &info);
549 dfeast_grci_(&ijob, &n,
reinterpret_cast<void*
>(&ze), zwork.GetDataVoid(),
550 workc.GetDataVoid(), zaq.GetDataVoid(), zbq.GetDataVoid(),
551 fpm.GetData(), &epsout, &loop, &Emid, &r,
552 &m0, lambda_cplx.GetDataVoid(), eigen_vec_cplx.GetDataVoid(),
553 &m, res.GetData(), &info);
556 void CallFeast(
int& ijob,
int n,
bool herm,
bool sym, complex<double>& ze,
557 Matrix<complex<double>, General, ColMajor>& work,
558 Matrix<complex<double>, General, ColMajor>& workc,
559 Vector<complex<double> >& aq, Vector<complex<double> >& bq,
560 Matrix<complex<double>, General, ColMajor>& zwork,
561 Vector<complex<double> >& zaq, Vector<complex<double> >& zbq,
562 IVect& fpm,
double& epsout,
int& loop,
563 double emin,
double emax, complex<double> Emid,
double r,
564 int m0, Vector<double>& lambda,
565 Matrix<complex<double>, General, ColMajor>& eigen_vectors,
566 Vector<complex<double> >& lambda_cplx,
567 Matrix<complex<double>, General, ColMajor>& eigen_vec_cplx,
568 int& m, Vector<double>& res,
int& info)
571 zfeast_hrci_(&ijob, &n,
reinterpret_cast<void*
>(&ze), work.GetDataVoid(),
572 workc.GetDataVoid(), aq.GetDataVoid(), bq.GetDataVoid(),
573 fpm.GetData(), &epsout, &loop, &emin, &emax, &m0,
574 lambda.GetData(), eigen_vectors.GetDataVoid(),
575 &m, res.GetData(), &info);
577 zfeast_srci_(&ijob, &n,
reinterpret_cast<void*
>(&ze), zwork.GetDataVoid(),
578 workc.GetDataVoid(), zaq.GetDataVoid(), zbq.GetDataVoid(),
579 fpm.GetData(), &epsout, &loop, &Emid, &r, &m0,
580 lambda_cplx.GetData(), eigen_vec_cplx.GetDataVoid(),
581 &m, res.GetData(), &info);
583 zfeast_grci_(&ijob, &n,
reinterpret_cast<void*
>(&ze), zwork.GetDataVoid(),
584 workc.GetDataVoid(), zaq.GetDataVoid(), zbq.GetDataVoid(),
585 fpm.GetData(), &epsout, &loop, &Emid, &r, &m0,
586 lambda_cplx.GetData(), eigen_vec_cplx.GetDataVoid(),
587 &m, res.GetData(), &info);
590 void CallFeastP(
int& ijob,
int p,
int n,
bool sym, complex<double>& ze,
591 Matrix<complex<double>, General, ColMajor>& zwork,
592 Matrix<complex<double>, General, ColMajor>& workc,
593 Vector<complex<double> >& zaq, Vector<complex<double> >& zbq,
594 IVect& fpm,
double& epsout,
int& loop, complex<double> Emid,
double r,
595 int m0, Vector<complex<double> >& lambda_cplx,
596 Matrix<complex<double>, General, ColMajor>& eigen_vec_cplx,
597 int& m, Vector<double>& res,
int& info)
600 zfeast_srcipev_(&ijob, &p, &n,
reinterpret_cast<void*
>(&ze),
601 zwork.GetDataVoid(), workc.GetDataVoid(),
602 zaq.GetDataVoid(), zbq.GetDataVoid(),
603 fpm.GetData(), &epsout, &loop, &Emid, &r, &m0,
604 lambda_cplx.GetData(), eigen_vec_cplx.GetDataVoid(),
605 &m, res.GetData(), &info);
607 zfeast_grcipev_(&ijob, &p, &n,
reinterpret_cast<void*
>(&ze),
608 zwork.GetDataVoid(), workc.GetDataVoid(),
609 zaq.GetDataVoid(), zbq.GetDataVoid(),
610 fpm.GetData(), &epsout, &loop, &Emid, &r, &m0,
611 lambda_cplx.GetData(), eigen_vec_cplx.GetDataVoid(),
612 &m, res.GetData(), &info);
615 #ifdef SELDON_WITH_VIRTUAL
617 void FindEigenvaluesFeast(PolynomialEigenProblem_Base<T>& var,
618 Vector<T>& eigen_values,
619 Vector<T>& lambda_imag,
620 Matrix<T, General, ColMajor>& eigen_vectors)
625 feastinit_(fpm.GetData());
628 if (var.GetPrintLevel() > 0)
631 double tol = var.GetStoppingCriterion();
632 fpm(2) = -log10(tol);
634 for (
int i = 64; i < fpm.GetM(); i++)
637 fpm(3) = var.GetNbMaximumIterations();
639 FeastParam& param = var.GetFeastParameters();
640 if (param.EstimateNumberEigenval())
643 #ifdef SELDON_WITH_MPI
644 MPI_Comm solver_comm;
645 MPI_Comm_dup(var.GetCommunicator(), &solver_comm);
647 fpm(8) = MPI_Comm_c2f(var.GetRciCommunicator());
648 fpm(48) = MPI_Comm_c2f(solver_comm);
649 fpm(58) = MPI_Comm_c2f(var.GetGlobalCommunicator());
654 int p = var.GetPolynomialDegree();
656 typedef typename ClassComplexType<T>::Tcplx Tcplx;
657 typedef typename ClassComplexType<T>::Treal Treal;
659 int m0 = var.GetNbAskedEigenvalues()+1;
661 Tcplx Emid = param.GetCircleCenterSpectrum();
662 Treal r = param.GetCircleRadiusSpectrum();
663 fpm(17) = param.GetRatioEllipseSpectrum();
664 fpm(18) = param.GetAngleEllipseSpectrum();
666 bool sym = var.IsSymmetricProblem();
672 if (param.GetTypeIntegration() >= 0)
673 fpm(15) = param.GetTypeIntegration();
675 if (param.GetNumOfQuadraturePoints() > 0)
676 fpm(7) = param.GetNumOfQuadraturePoints();
681 Matrix<Tcplx, General, ColMajor> zwork;
682 Matrix<Tcplx, General, ColMajor> workc;
683 Vector<Tcplx> xc(n), yc(n);
684 Vector<Tcplx> zaq, zbq;
686 workc.Reallocate(n, m0);
687 zwork.Reallocate(n, m0);
688 zaq.Reallocate(p*p*m0*m0);
689 zbq.Reallocate(p*p*m0*m0);
692 workc.Zero(); xc.Zero(); yc.Zero();
693 zaq.Zero(); zbq.Zero();
696 Treal epsout(0); Tcplx ze;
697 Vector<Treal> res(m0);
700 Matrix<Tcplx, General, ColMajor> eigen_vec_cplx;
701 Vector<Treal> lambda; Vector<Tcplx> lambda_cplx;
703 eigen_vec_cplx.Reallocate(n, m0);
704 lambda_cplx.Reallocate(m0);
707 eigen_vectors.Zero(); lambda.Zero();
708 lambda_cplx.Zero(); eigen_vec_cplx.Zero();
710 Vector<Tcplx> coef(p+1); coef.Zero();
713 int ijob = -1, info = 0;
717 CallFeastP(ijob, p, n, sym, ze, zwork, workc, zaq, zbq,
718 fpm, epsout, loop, Emid, r, m0, lambda_cplx, eigen_vec_cplx,
724 SetComplexOne(coef(0));
725 for (
int i = 0; i < p; i++)
726 coef(i+1) = coef(i)*ze;
728 var.FactorizeOperator(coef);
733 for (
int k = 0; k < fpm(22); k++)
735 for (
int i = 0; i < n; i++)
738 var.SolveOperator(SeldonNoTrans, xc, yc);
739 var.IncrementLinearSolves();
741 for (
int i = 0; i < n; i++)
753 for (
int k = 0; k < fpm(22); k++)
755 for (
int i = 0; i < n; i++)
759 var.SolveOperator(SeldonTrans, xc, yc);
760 var.IncrementLinearSolves();
763 for (
int i = 0; i < n; i++)
770 int i = fpm(23), j = fpm(23) + fpm(24)-1;
771 for (
int k = i; k <= j; k++)
773 GetCol(eigen_vec_cplx, k-1, xc);
775 var.MltOperator(fpm(56)-1, SeldonNoTrans, xc, yc);
776 var.IncrementProdMatVect();
778 SetCol(yc, k-1, zwork);
784 int i = fpm(33), j = fpm(33) + fpm(34)-1;
785 for (
int k = i; k <= j; k++)
787 GetCol(eigen_vec_cplx, k-1, xc);
789 var.MltOperator(fpm(56)-1, SeldonConjTrans, xc, yc);
790 var.IncrementProdMatVect();
792 SetCol(yc, k-1, zwork);
796 zwork.Clear(); workc.Clear(); zaq.Clear(); zbq.Clear();
798 var.FillComplexEigenvectors(m, Emid, tol, lambda_cplx, eigen_vec_cplx,
799 eigen_values, lambda_imag, eigen_vectors);
801 eigen_vec_cplx.Clear(); lambda_cplx.Clear();
803 var.DistributeEigenvectors(eigen_vectors);