Actual source code: sinvert.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:       Implements the shift-and-invert technique for eigenvalue problems.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/stimpl.h>

 28: PetscErrorCode STApply_Sinvert(ST st,Vec x,Vec y)
 29: {

 33:   if (st->nmat>1) {
 34:     /* generalized eigenproblem: y = (A - sB)^-1 B x */
 35:     MatMult(st->T[0],x,st->w);
 36:     STMatSolve(st,st->w,y);
 37:   } else {
 38:     /* standard eigenproblem: y = (A - sI)^-1 x */
 39:     STMatSolve(st,x,y);
 40:   }
 41:   return(0);
 42: }

 46: PetscErrorCode STApplyTranspose_Sinvert(ST st,Vec x,Vec y)
 47: {

 51:   if (st->nmat>1) {
 52:     /* generalized eigenproblem: y = B^T (A - sB)^-T x */
 53:     STMatSolveTranspose(st,x,st->w);
 54:     MatMultTranspose(st->T[0],st->w,y);
 55:   } else {
 56:     /* standard eigenproblem: y = (A - sI)^-T x */
 57:     STMatSolveTranspose(st,x,y);
 58:   }
 59:   return(0);
 60: }

 64: PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 65: {
 66:   PetscInt    j;
 67: #if !defined(PETSC_USE_COMPLEX)
 68:   PetscScalar t;
 69: #endif

 72: #if !defined(PETSC_USE_COMPLEX)
 73:   for (j=0;j<n;j++) {
 74:     if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
 75:     else {
 76:       t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
 77:       eigr[j] = eigr[j] / t + st->sigma;
 78:       eigi[j] = - eigi[j] / t;
 79:     }
 80:   }
 81: #else
 82:   for (j=0;j<n;j++) {
 83:     eigr[j] = 1.0 / eigr[j] + st->sigma;
 84:   }
 85: #endif
 86:   return(0);
 87: }

 91: PetscErrorCode STPostSolve_Sinvert(ST st)
 92: {

 96:   if (st->shift_matrix == ST_MATMODE_INPLACE) {
 97:     if (st->nmat>1) {
 98:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
 99:     } else {
100:       MatShift(st->A[0],st->sigma);
101:     }
102:     st->Astate[0] = ((PetscObject)st->A[0])->state;
103:     st->state = ST_STATE_INITIAL;
104:   }
105:   return(0);
106: }

110: PetscErrorCode STSetUp_Sinvert(ST st)
111: {
113:   PetscInt       k,nc,nmat=PetscMax(st->nmat,2);
114:   PetscScalar    *coeffs=NULL;

117:   if (st->nmat>1) {
118:     ST_AllocateWorkVec(st);
119:   }
120:   /* if the user did not set the shift, use the target value */
121:   if (!st->sigma_set) st->sigma = st->defsigma;
122:   if (st->transform) {
123:     if (nmat>2) {
124:       nc = (nmat*(nmat+1))/2;
125:       PetscMalloc(nc*sizeof(PetscScalar),&coeffs);
126:       /* Compute coeffs */
127:       STCoeffs_Monomial(st,coeffs);
128:     }
129:     /* T[0] = A_n */
130:     k = nmat-1;
131:     PetscObjectReference((PetscObject)st->A[k]);
132:     MatDestroy(&st->T[0]);
133:     st->T[0] = st->A[k];
134:     for (k=1;k<nmat;k++) {
135:       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]);
136:     }
137:     if (nmat>2) { PetscFree(coeffs); }
138:     PetscObjectReference((PetscObject)st->T[nmat-1]);
139:     MatDestroy(&st->P);
140:     st->P = st->T[nmat-1];
141:   } else {
142:     for (k=0;k<nmat;k++) {
143:       PetscObjectReference((PetscObject)st->A[k]);
144:       MatDestroy(&st->T[k]);
145:       st->T[k] = st->A[k];
146:     }
147:   } 
148:   if (st->P) {
149:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
150:     STCheckFactorPackage(st);
151:     KSPSetOperators(st->ksp,st->P,st->P);
152:     KSPSetErrorIfNotConverged(st->ksp,PETSC_TRUE);
153:     KSPSetUp(st->ksp);
154:   }
155:   return(0);
156: }

160: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
161: {
163:   PetscInt       nmat=PetscMax(st->nmat,2),k,nc;
164:   PetscScalar    *coeffs=NULL;

167:   if (st->transform) {
168:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
169:       nc = (nmat*(nmat+1))/2;
170:       PetscMalloc(nc*sizeof(PetscScalar),&coeffs);
171:       /* Compute coeffs */
172:       STCoeffs_Monomial(st,coeffs);
173:     }
174:     for (k=1;k<nmat;k++) {
175:       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]);
176:     }
177:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
178:       PetscFree(coeffs);
179:     }
180:     if (st->P!=st->T[nmat-1]) {
181:       MatDestroy(&st->P);
182:       st->P = st->T[nmat-1];
183:       PetscObjectReference((PetscObject)st->P);
184:     }
185:   }
186:   if (st->P) {
187:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
188:     KSPSetOperators(st->ksp,st->P,st->P);
189:     KSPSetUp(st->ksp);
190:   }
191:   return(0);
192: }

196: PetscErrorCode STSetFromOptions_Sinvert(PetscOptionItems *PetscOptionsObject,ST st)
197: {
199:   PC             pc;
200:   PCType         pctype;
201:   KSPType        ksptype;

204:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
205:   KSPGetPC(st->ksp,&pc);
206:   KSPGetType(st->ksp,&ksptype);
207:   PCGetType(pc,&pctype);
208:   if (!pctype && !ksptype) {
209:     if (st->shift_matrix == ST_MATMODE_SHELL) {
210:       /* in shell mode use GMRES with Jacobi as the default */
211:       KSPSetType(st->ksp,KSPGMRES);
212:       PCSetType(pc,PCJACOBI);
213:     } else {
214:       /* use direct solver as default */
215:       KSPSetType(st->ksp,KSPPREONLY);
216:       PCSetType(pc,PCLU);
217:     }
218:   }
219:   return(0);
220: }

224: PETSC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
225: {
227:   st->ops->apply           = STApply_Sinvert;
228:   st->ops->getbilinearform = STGetBilinearForm_Default;
229:   st->ops->applytrans      = STApplyTranspose_Sinvert;
230:   st->ops->postsolve       = STPostSolve_Sinvert;
231:   st->ops->backtransform   = STBackTransform_Sinvert;
232:   st->ops->setup           = STSetUp_Sinvert;
233:   st->ops->setshift        = STSetShift_Sinvert;
234:   st->ops->setfromoptions  = STSetFromOptions_Sinvert;
235:   st->ops->checknullspace  = STCheckNullSpace_Default;
236:   return(0);
237: }