Actual source code: shift.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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-2014, 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->setupcalled = 0;
 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:     st->T[k] = st->A[k];
115:     for (k=0;k<nmat-1;k++) {
116:       STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PETSC_TRUE,&st->T[k]);
117:     }
118:      if (nmat>2) { PetscFree(coeffs); }
119:   } else {
120:     for (k=0;k<nmat;k++) {
121:       PetscObjectReference((PetscObject)st->A[k]);
122:       st->T[k] = st->A[k];
123:     }
124:   }
125:   if (nmat==2 || (nmat>2 && st->transform)) {
126:     st->P = st->T[nmat-1];
127:     PetscObjectReference((PetscObject)st->P);
128:   }
129:   if (st->P) {
130:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
131:     KSPSetOperators(st->ksp,st->P,st->P);
132:     KSPSetUp(st->ksp);
133:   }
134:   return(0);
135: }

139: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
140: {
142:   PetscInt       k,nc,nmat=PetscMax(st->nmat,2);
143:   PetscScalar    *coeffs=NULL;

146:   /* Nothing to be done if STSetUp has not been called yet */
147:   if (!st->setupcalled) return(0);

149:   if (nmat<3 || st->transform) {
150:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
151:       nc = (nmat*(nmat+1))/2;
152:       PetscMalloc(nc*sizeof(PetscScalar),&coeffs);
153:       /* Compute coeffs */
154:       STCoeffs_Monomial(st,coeffs);
155:     }
156:     for (k=0;k<nmat-1;k++) {
157:       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]);
158:     }
159:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
160:         PetscFree(coeffs);
161:     }
162:   }
163:   return(0);
164: }

168: PetscErrorCode STSetFromOptions_Shift(ST st)
169: {
171:   PC             pc;
172:   PCType         pctype;
173:   KSPType        ksptype;

176:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
177:   KSPGetPC(st->ksp,&pc);
178:   KSPGetType(st->ksp,&ksptype);
179:   PCGetType(pc,&pctype);
180:   if (!pctype && !ksptype) {
181:     if (st->shift_matrix == ST_MATMODE_SHELL) {
182:       /* in shell mode use GMRES with Jacobi as the default */
183:       KSPSetType(st->ksp,KSPGMRES);
184:       PCSetType(pc,PCJACOBI);
185:     } else {
186:       /* use direct solver as default */
187:       KSPSetType(st->ksp,KSPPREONLY);
188:       PCSetType(pc,PCREDUNDANT);
189:     }
190:   }
191:   return(0);
192: }

196: PETSC_EXTERN PetscErrorCode STCreate_Shift(ST st)
197: {
199:   st->ops->apply           = STApply_Shift;
200:   st->ops->getbilinearform = STGetBilinearForm_Default;
201:   st->ops->applytrans      = STApplyTranspose_Shift;
202:   st->ops->postsolve       = STPostSolve_Shift;
203:   st->ops->backtransform   = STBackTransform_Shift;
204:   st->ops->setfromoptions  = STSetFromOptions_Shift;
205:   st->ops->setup           = STSetUp_Shift;
206:   st->ops->setshift        = STSetShift_Shift;
207:   return(0);
208: }