Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 : /*
11 : This file implements a wrapper to the LAPACK eigenvalue subroutines.
12 : Generalized problems are transformed to standard ones only if necessary.
13 : */
14 :
15 : #include <slepc/private/epsimpl.h>
16 :
17 28 : static PetscErrorCode EPSSetUp_LAPACK(EPS eps)
18 : {
19 28 : int ierra,ierrb;
20 28 : PetscBool isshift,flg,denseok=PETSC_FALSE;
21 28 : Mat A,B,OP,shell,Ar,Br,Adense=NULL,Bdense=NULL,Ads,Bds;
22 28 : PetscScalar shift;
23 28 : PetscInt nmat;
24 28 : KSP ksp;
25 28 : PC pc;
26 :
27 28 : PetscFunctionBegin;
28 28 : eps->ncv = eps->n;
29 28 : if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
30 28 : if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
31 28 : if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
32 28 : PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
33 28 : EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
34 28 : EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
35 28 : PetscCall(EPSAllocateSolution(eps,0));
36 :
37 : /* attempt to get dense representations of A and B separately */
38 28 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
39 28 : if (isshift) {
40 27 : PetscCall(STGetNumMatrices(eps->st,&nmat));
41 27 : PetscCall(STGetMatrix(eps->st,0,&A));
42 27 : PetscCall(MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg));
43 27 : if (flg) {
44 24 : PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler,NULL));
45 24 : ierra = MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
46 24 : if (!ierra) ierra |= MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense);
47 24 : ierra |= MatDestroy(&Ar);
48 24 : PetscCall(PetscPopErrorHandler());
49 : } else ierra = 1;
50 27 : if (nmat>1) {
51 1 : PetscCall(STGetMatrix(eps->st,1,&B));
52 1 : PetscCall(MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg));
53 1 : if (flg) {
54 1 : PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler,NULL));
55 1 : ierrb = MatCreateRedundantMatrix(B,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Br);
56 1 : if (!ierrb) ierrb |= MatConvert(Br,MATSEQDENSE,MAT_INITIAL_MATRIX,&Bdense);
57 1 : ierrb |= MatDestroy(&Br);
58 1 : PetscCall(PetscPopErrorHandler());
59 : } else ierrb = 1;
60 : } else ierrb = 0;
61 27 : denseok = PetscNot(ierra || ierrb);
62 : }
63 :
64 : /* setup DS */
65 27 : if (denseok) {
66 24 : if (eps->isgeneralized) {
67 1 : if (eps->ishermitian) {
68 1 : if (eps->ispositive) PetscCall(DSSetType(eps->ds,DSGHEP));
69 0 : else PetscCall(DSSetType(eps->ds,DSGNHEP)); /* TODO: should be DSGHIEP */
70 0 : } else PetscCall(DSSetType(eps->ds,DSGNHEP));
71 : } else {
72 23 : if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
73 1 : else PetscCall(DSSetType(eps->ds,DSNHEP));
74 : }
75 4 : } else PetscCall(DSSetType(eps->ds,DSNHEP));
76 28 : PetscCall(DSAllocate(eps->ds,eps->ncv));
77 28 : PetscCall(DSSetDimensions(eps->ds,eps->ncv,0,0));
78 :
79 28 : if (denseok) {
80 24 : PetscCall(STGetShift(eps->st,&shift));
81 24 : if (shift != 0.0) {
82 0 : if (nmat>1) PetscCall(MatAXPY(Adense,-shift,Bdense,SAME_NONZERO_PATTERN));
83 0 : else PetscCall(MatShift(Adense,-shift));
84 : }
85 : /* use dummy pc and ksp to avoid problems when B is not positive definite */
86 24 : PetscCall(STGetKSP(eps->st,&ksp));
87 24 : PetscCall(KSPSetType(ksp,KSPPREONLY));
88 24 : PetscCall(KSPGetPC(ksp,&pc));
89 24 : PetscCall(PCSetType(pc,PCNONE));
90 : } else {
91 4 : PetscCall(PetscInfo(eps,"Using slow explicit operator\n"));
92 4 : PetscCall(STGetOperator(eps->st,&shell));
93 4 : PetscCall(MatComputeOperator(shell,MATDENSE,&OP));
94 4 : PetscCall(STRestoreOperator(eps->st,&shell));
95 4 : PetscCall(MatDestroy(&Adense));
96 4 : PetscCall(MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Adense));
97 4 : PetscCall(MatDestroy(&OP));
98 : }
99 :
100 : /* fill DS matrices */
101 28 : PetscCall(DSGetMat(eps->ds,DS_MAT_A,&Ads));
102 28 : PetscCall(MatCopy(Adense,Ads,SAME_NONZERO_PATTERN));
103 28 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&Ads));
104 28 : if (denseok && eps->isgeneralized) {
105 1 : PetscCall(DSGetMat(eps->ds,DS_MAT_B,&Bds));
106 1 : PetscCall(MatCopy(Bdense,Bds,SAME_NONZERO_PATTERN));
107 1 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&Bds));
108 : }
109 28 : PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
110 28 : PetscCall(MatDestroy(&Adense));
111 28 : PetscCall(MatDestroy(&Bdense));
112 28 : PetscFunctionReturn(PETSC_SUCCESS);
113 : }
114 :
115 28 : static PetscErrorCode EPSSolve_LAPACK(EPS eps)
116 : {
117 28 : PetscInt n=eps->n,i,low,high;
118 28 : PetscScalar *array,*pX,*pY;
119 28 : Vec v,w;
120 :
121 28 : PetscFunctionBegin;
122 28 : PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
123 28 : PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
124 28 : PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
125 :
126 : /* right eigenvectors */
127 28 : PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
128 28 : PetscCall(DSGetArray(eps->ds,DS_MAT_X,&pX));
129 3332 : for (i=0;i<eps->ncv;i++) {
130 3304 : PetscCall(BVGetColumn(eps->V,i,&v));
131 3304 : PetscCall(VecGetOwnershipRange(v,&low,&high));
132 3304 : PetscCall(VecGetArray(v,&array));
133 3304 : PetscCall(PetscArraycpy(array,pX+i*n+low,high-low));
134 3304 : PetscCall(VecRestoreArray(v,&array));
135 3304 : PetscCall(BVRestoreColumn(eps->V,i,&v));
136 : }
137 28 : PetscCall(DSRestoreArray(eps->ds,DS_MAT_X,&pX));
138 :
139 : /* left eigenvectors */
140 28 : if (eps->twosided) {
141 1 : PetscCall(DSVectors(eps->ds,DS_MAT_Y,NULL,NULL));
142 1 : PetscCall(DSGetArray(eps->ds,DS_MAT_Y,&pY));
143 101 : for (i=0;i<eps->ncv;i++) {
144 100 : PetscCall(BVGetColumn(eps->W,i,&w));
145 100 : PetscCall(VecGetOwnershipRange(w,&low,&high));
146 100 : PetscCall(VecGetArray(w,&array));
147 100 : PetscCall(PetscArraycpy(array,pY+i*n+low,high-low));
148 100 : PetscCall(VecRestoreArray(w,&array));
149 100 : PetscCall(BVRestoreColumn(eps->W,i,&w));
150 : }
151 1 : PetscCall(DSRestoreArray(eps->ds,DS_MAT_Y,&pY));
152 : }
153 :
154 28 : eps->nconv = eps->ncv;
155 28 : eps->its = 1;
156 28 : eps->reason = EPS_CONVERGED_TOL;
157 28 : PetscFunctionReturn(PETSC_SUCCESS);
158 : }
159 :
160 17 : SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
161 : {
162 17 : PetscFunctionBegin;
163 17 : eps->useds = PETSC_TRUE;
164 17 : eps->categ = EPS_CATEGORY_OTHER;
165 :
166 17 : eps->ops->solve = EPSSolve_LAPACK;
167 17 : eps->ops->setup = EPSSetUp_LAPACK;
168 17 : eps->ops->setupsort = EPSSetUpSort_Default;
169 17 : eps->ops->backtransform = EPSBackTransform_Default;
170 17 : PetscFunctionReturn(PETSC_SUCCESS);
171 : }
|