Actual source code: sinvert.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Implements the shift-and-invert technique for eigenvalue problems
12: */
14: #include <slepc/private/stimpl.h>
16: PetscErrorCode STApply_Sinvert(ST st,Vec x,Vec y)
17: {
21: if (st->nmat>1) {
22: /* generalized eigenproblem: y = (A - sB)^-1 B x */
23: MatMult(st->T[0],x,st->work[0]);
24: STMatSolve(st,st->work[0],y);
25: } else {
26: /* standard eigenproblem: y = (A - sI)^-1 x */
27: STMatSolve(st,x,y);
28: }
29: return(0);
30: }
32: PetscErrorCode STApplyTranspose_Sinvert(ST st,Vec x,Vec y)
33: {
37: if (st->nmat>1) {
38: /* generalized eigenproblem: y = B^T (A - sB)^-T x */
39: STMatSolveTranspose(st,x,st->work[0]);
40: MatMultTranspose(st->T[0],st->work[0],y);
41: } else {
42: /* standard eigenproblem: y = (A - sI)^-T x */
43: STMatSolveTranspose(st,x,y);
44: }
45: return(0);
46: }
48: PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
49: {
50: PetscInt j;
51: #if !defined(PETSC_USE_COMPLEX)
52: PetscScalar t;
53: #endif
56: #if !defined(PETSC_USE_COMPLEX)
57: for (j=0;j<n;j++) {
58: if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
59: else {
60: t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
61: eigr[j] = eigr[j] / t + st->sigma;
62: eigi[j] = - eigi[j] / t;
63: }
64: }
65: #else
66: for (j=0;j<n;j++) {
67: eigr[j] = 1.0 / eigr[j] + st->sigma;
68: }
69: #endif
70: return(0);
71: }
73: PetscErrorCode STPostSolve_Sinvert(ST st)
74: {
78: if (st->shift_matrix == ST_MATMODE_INPLACE) {
79: if (st->nmat>1) {
80: MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
81: } else {
82: MatShift(st->A[0],st->sigma);
83: }
84: st->Astate[0] = ((PetscObject)st->A[0])->state;
85: st->state = ST_STATE_INITIAL;
86: }
87: return(0);
88: }
90: PetscErrorCode STSetUp_Sinvert(ST st)
91: {
93: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
94: PetscScalar *coeffs=NULL;
97: if (st->nmat>1) {
98: STSetWorkVecs(st,1);
99: }
100: /* if the user did not set the shift, use the target value */
101: if (!st->sigma_set) st->sigma = st->defsigma;
102: if (st->transform) {
103: if (nmat>2) {
104: nc = (nmat*(nmat+1))/2;
105: PetscMalloc1(nc,&coeffs);
106: /* Compute coeffs */
107: STCoeffs_Monomial(st,coeffs);
108: }
109: /* T[0] = A_n */
110: k = nmat-1;
111: PetscObjectReference((PetscObject)st->A[k]);
112: MatDestroy(&st->T[0]);
113: st->T[0] = st->A[k];
114: for (k=1;k<nmat;k++) {
115: STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,nmat-k-1,coeffs?coeffs+(k*(k+1))/2:NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[k]);
116: }
117: if (nmat>2) { PetscFree(coeffs); }
118: PetscObjectReference((PetscObject)st->T[nmat-1]);
119: MatDestroy(&st->P);
120: st->P = st->T[nmat-1];
121: } else {
122: for (k=0;k<nmat;k++) {
123: PetscObjectReference((PetscObject)st->A[k]);
124: MatDestroy(&st->T[k]);
125: st->T[k] = st->A[k];
126: }
127: }
128: if (st->P) {
129: if (!st->ksp) { STGetKSP(st,&st->ksp); }
130: STCheckFactorPackage(st);
131: KSPSetOperators(st->ksp,st->P,st->P);
132: KSPSetUp(st->ksp);
133: }
134: return(0);
135: }
137: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
138: {
140: PetscInt nmat=PetscMax(st->nmat,2),k,nc;
141: PetscScalar *coeffs=NULL;
144: if (st->transform) {
145: if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
146: nc = (nmat*(nmat+1))/2;
147: PetscMalloc1(nc,&coeffs);
148: /* Compute coeffs */
149: STCoeffs_Monomial(st,coeffs);
150: }
151: for (k=1;k<nmat;k++) {
152: STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,nmat-k-1,coeffs?coeffs+(k*(k+1))/2:NULL,PETSC_FALSE,&st->T[k]);
153: }
154: if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
155: PetscFree(coeffs);
156: }
157: if (st->P!=st->T[nmat-1]) {
158: PetscObjectReference((PetscObject)st->T[nmat-1]);
159: MatDestroy(&st->P);
160: st->P = st->T[nmat-1];
161: }
162: }
163: if (st->P) {
164: if (!st->ksp) { STGetKSP(st,&st->ksp); }
165: KSPSetOperators(st->ksp,st->P,st->P);
166: KSPSetUp(st->ksp);
167: }
168: return(0);
169: }
171: SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
172: {
174: st->ops->apply = STApply_Sinvert;
175: st->ops->getbilinearform = STGetBilinearForm_Default;
176: st->ops->applytrans = STApplyTranspose_Sinvert;
177: st->ops->postsolve = STPostSolve_Sinvert;
178: st->ops->backtransform = STBackTransform_Sinvert;
179: st->ops->setup = STSetUp_Sinvert;
180: st->ops->setshift = STSetShift_Sinvert;
181: st->ops->checknullspace = STCheckNullSpace_Default;
182: st->ops->setdefaultksp = STSetDefaultKSP_Default;
183: return(0);
184: }