Actual source code: sinvert.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 17: {
 18:   PetscInt    j;
 19: #if !defined(PETSC_USE_COMPLEX)
 20:   PetscScalar t;
 21: #endif

 24: #if !defined(PETSC_USE_COMPLEX)
 25:   for (j=0;j<n;j++) {
 26:     if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
 27:     else {
 28:       t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
 29:       eigr[j] = eigr[j] / t + st->sigma;
 30:       eigi[j] = - eigi[j] / t;
 31:     }
 32:   }
 33: #else
 34:   for (j=0;j<n;j++) {
 35:     eigr[j] = 1.0 / eigr[j] + st->sigma;
 36:   }
 37: #endif
 38:   return(0);
 39: }

 41: PetscErrorCode STPostSolve_Sinvert(ST st)
 42: {

 46:   if (st->matmode == ST_MATMODE_INPLACE) {
 47:     if (st->nmat>1) {
 48:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
 49:     } else {
 50:       MatShift(st->A[0],st->sigma);
 51:     }
 52:     st->Astate[0] = ((PetscObject)st->A[0])->state;
 53:     st->state   = ST_STATE_INITIAL;
 54:     st->opready = PETSC_FALSE;
 55:   }
 56:   return(0);
 57: }

 59: /*
 60:    Operator (sinvert):
 61:                Op               P         M
 62:    if nmat=1:  (A-sI)^-1        A-sI      NULL
 63:    if nmat=2:  (A-sB)^-1 B      A-sB      B
 64: */
 65: PetscErrorCode STComputeOperator_Sinvert(ST st)
 66: {

 70:   /* if the user did not set the shift, use the target value */
 71:   if (!st->sigma_set) st->sigma = st->defsigma;
 72:   PetscObjectReference((PetscObject)st->A[1]);
 73:   MatDestroy(&st->T[0]);
 74:   st->T[0] = st->A[1];
 75:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[1]);
 76:   PetscObjectReference((PetscObject)st->T[1]);
 77:   MatDestroy(&st->P);
 78:   st->P = st->T[1];
 79:   st->M = (st->nmat>1)? st->T[0]: NULL;
 80:   return(0);
 81: }

 83: PetscErrorCode STSetUp_Sinvert(ST st)
 84: {
 86:   PetscInt       k,nc,nmat=st->nmat;
 87:   PetscScalar    *coeffs=NULL;

 90:   if (nmat>1) {
 91:     STSetWorkVecs(st,1);
 92:   }
 93:   /* if the user did not set the shift, use the target value */
 94:   if (!st->sigma_set) st->sigma = st->defsigma;
 95:   if (nmat>2) {  /* set-up matrices for polynomial eigenproblems */
 96:     if (st->transform) {
 97:       nc = (nmat*(nmat+1))/2;
 98:       PetscMalloc1(nc,&coeffs);
 99:       /* Compute coeffs */
100:       STCoeffs_Monomial(st,coeffs);
101:       /* T[0] = A_n */
102:       k = nmat-1;
103:       PetscObjectReference((PetscObject)st->A[k]);
104:       MatDestroy(&st->T[0]);
105:       st->T[0] = st->A[k];
106:       for (k=1;k<nmat;k++) {
107:         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]);
108:       }
109:       PetscFree(coeffs);
110:       PetscObjectReference((PetscObject)st->T[nmat-1]);
111:       MatDestroy(&st->P);
112:       st->P = st->T[nmat-1];
113:       STKSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
114:     } else {
115:       for (k=0;k<nmat;k++) {
116:         PetscObjectReference((PetscObject)st->A[k]);
117:         MatDestroy(&st->T[k]);
118:         st->T[k] = st->A[k];
119:       }
120:     }
121:   }
122:   if (st->P) {
123:     KSPSetUp(st->ksp);
124:   }
125:   return(0);
126: }

128: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
129: {
131:   PetscInt       nmat=PetscMax(st->nmat,2),k,nc;
132:   PetscScalar    *coeffs=NULL;

135:   if (st->transform) {
136:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
137:       nc = (nmat*(nmat+1))/2;
138:       PetscMalloc1(nc,&coeffs);
139:       /* Compute coeffs */
140:       STCoeffs_Monomial(st,coeffs);
141:     }
142:     for (k=1;k<nmat;k++) {
143:       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]);
144:     }
145:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
146:       PetscFree(coeffs);
147:     }
148:     if (st->P!=st->T[nmat-1]) {
149:       PetscObjectReference((PetscObject)st->T[nmat-1]);
150:       MatDestroy(&st->P);
151:       st->P = st->T[nmat-1];
152:     }
153:   }
154:   if (st->P) {
155:     STKSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
156:     KSPSetUp(st->ksp);
157:   }
158:   return(0);
159: }

161: SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
162: {
164:   st->usesksp = PETSC_TRUE;

166:   st->ops->apply           = STApply_Generic;
167:   st->ops->applytrans      = STApplyTranspose_Generic;
168:   st->ops->backtransform   = STBackTransform_Sinvert;
169:   st->ops->setshift        = STSetShift_Sinvert;
170:   st->ops->getbilinearform = STGetBilinearForm_Default;
171:   st->ops->setup           = STSetUp_Sinvert;
172:   st->ops->computeoperator = STComputeOperator_Sinvert;
173:   st->ops->postsolve       = STPostSolve_Sinvert;
174:   st->ops->checknullspace  = STCheckNullSpace_Default;
175:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
176:   return(0);
177: }