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