Actual source code: shift.c
slepc-3.7.0 2016-05-16
1: /*
2: Shift spectral transformation, applies (A + sigma I) as operator, or
3: inv(B)(A + sigma B) for generalized problems
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <slepc/private/stimpl.h>
29: PetscErrorCode STApply_Shift(ST st,Vec x,Vec y)
30: {
34: if (st->nmat>1) {
35: /* generalized eigenproblem: y = B^-1 (A - sB) x */
36: MatMult(st->T[0],x,st->w);
37: STMatSolve(st,st->w,y);
38: } else {
39: /* standard eigenproblem: y = (A - sI) x */
40: MatMult(st->T[0],x,y);
41: }
42: return(0);
43: }
47: PetscErrorCode STApplyTranspose_Shift(ST st,Vec x,Vec y)
48: {
52: if (st->nmat>1) {
53: /* generalized eigenproblem: y = (A - sB)^T B^-T x */
54: STMatSolveTranspose(st,x,st->w);
55: MatMultTranspose(st->T[0],st->w,y);
56: } else {
57: /* standard eigenproblem: y = (A^T - sI) x */
58: MatMultTranspose(st->T[0],x,y);
59: }
60: return(0);
61: }
65: PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
66: {
67: PetscInt j;
70: for (j=0;j<n;j++) {
71: eigr[j] += st->sigma;
72: }
73: return(0);
74: }
78: PetscErrorCode STPostSolve_Shift(ST st)
79: {
83: if (st->shift_matrix == ST_MATMODE_INPLACE) {
84: if (st->nmat>1) {
85: MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
86: } else {
87: MatShift(st->A[0],st->sigma);
88: }
89: st->Astate[0] = ((PetscObject)st->A[0])->state;
90: st->state = ST_STATE_INITIAL;
91: }
92: return(0);
93: }
97: PetscErrorCode STSetUp_Shift(ST st)
98: {
100: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
101: PetscScalar *coeffs=NULL;
104: if (nmat<3 || st->transform) {
105: if (nmat>2) {
106: nc = (nmat*(nmat+1))/2;
107: PetscMalloc(nc*sizeof(PetscScalar),&coeffs);
108: /* Compute coeffs */
109: STCoeffs_Monomial(st,coeffs);
110: }
111: /* T[n] = A_n */
112: k = nmat-1;
113: PetscObjectReference((PetscObject)st->A[k]);
114: MatDestroy(&st->T[k]);
115: st->T[k] = st->A[k];
116: for (k=0;k<nmat-1;k++) {
117: STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[k]);
118: }
119: if (nmat>2) { PetscFree(coeffs); }
120: } else {
121: for (k=0;k<nmat;k++) {
122: PetscObjectReference((PetscObject)st->A[k]);
123: MatDestroy(&st->T[k]);
124: st->T[k] = st->A[k];
125: }
126: }
127: if (nmat>=2 && st->transform) {
128: PetscObjectReference((PetscObject)st->T[nmat-1]);
129: MatDestroy(&st->P);
130: st->P = st->T[nmat-1];
131: }
132: if (st->P) {
133: if (!st->ksp) { STGetKSP(st,&st->ksp); }
134: STCheckFactorPackage(st);
135: KSPSetOperators(st->ksp,st->P,st->P);
136: KSPSetErrorIfNotConverged(st->ksp,PETSC_TRUE);
137: KSPSetUp(st->ksp);
138: }
139: return(0);
140: }
144: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
145: {
147: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
148: PetscScalar *coeffs=NULL;
151: if (st->transform) {
152: if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
153: nc = (nmat*(nmat+1))/2;
154: PetscMalloc(nc*sizeof(PetscScalar),&coeffs);
155: /* Compute coeffs */
156: STCoeffs_Monomial(st,coeffs);
157: }
158: for (k=0;k<nmat-1;k++) {
159: STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PETSC_FALSE,&st->T[k]);
160: }
161: if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
162: PetscFree(coeffs);
163: }
164: }
165: return(0);
166: }
170: PetscErrorCode STSetFromOptions_Shift(PetscOptionItems *PetscOptionsObject,ST st)
171: {
173: PC pc;
174: PCType pctype;
175: KSPType ksptype;
178: if (!st->ksp) { STGetKSP(st,&st->ksp); }
179: KSPGetPC(st->ksp,&pc);
180: KSPGetType(st->ksp,&ksptype);
181: PCGetType(pc,&pctype);
182: if (!pctype && !ksptype) {
183: if (st->shift_matrix == ST_MATMODE_SHELL) {
184: /* in shell mode use GMRES with Jacobi as the default */
185: KSPSetType(st->ksp,KSPGMRES);
186: PCSetType(pc,PCJACOBI);
187: } else {
188: /* use direct solver as default */
189: KSPSetType(st->ksp,KSPPREONLY);
190: PCSetType(pc,PCLU);
191: }
192: }
193: return(0);
194: }
198: PETSC_EXTERN PetscErrorCode STCreate_Shift(ST st)
199: {
201: st->ops->apply = STApply_Shift;
202: st->ops->getbilinearform = STGetBilinearForm_Default;
203: st->ops->applytrans = STApplyTranspose_Shift;
204: st->ops->postsolve = STPostSolve_Shift;
205: st->ops->backtransform = STBackTransform_Shift;
206: st->ops->setfromoptions = STSetFromOptions_Shift;
207: st->ops->setup = STSetUp_Shift;
208: st->ops->setshift = STSetShift_Shift;
209: return(0);
210: }