Actual source code: sinvert.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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

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

 40: PetscErrorCode STPostSolve_Sinvert(ST st)
 41: {
 42:   if (st->matmode == ST_MATMODE_INPLACE) {
 43:     if (st->nmat>1) MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
 44:     else MatShift(st->A[0],st->sigma);
 45:     st->Astate[0] = ((PetscObject)st->A[0])->state;
 46:     st->state   = ST_STATE_INITIAL;
 47:     st->opready = PETSC_FALSE;
 48:   }
 49:   PetscFunctionReturn(0);
 50: }

 52: /*
 53:    Operator (sinvert):
 54:                Op               P         M
 55:    if nmat=1:  (A-sI)^-1        A-sI      NULL
 56:    if nmat=2:  (A-sB)^-1 B      A-sB      B
 57: */
 58: PetscErrorCode STComputeOperator_Sinvert(ST st)
 59: {
 60:   /* if the user did not set the shift, use the target value */
 61:   if (!st->sigma_set) st->sigma = st->defsigma;
 62:   PetscObjectReference((PetscObject)st->A[1]);
 63:   MatDestroy(&st->T[0]);
 64:   st->T[0] = st->A[1];
 65:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]);
 66:   PetscObjectReference((PetscObject)st->T[1]);
 67:   MatDestroy(&st->P);
 68:   st->P = st->T[1];
 69:   st->M = (st->nmat>1)? st->T[0]: NULL;
 70:   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
 71:     STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat);
 72:   }
 73:   PetscFunctionReturn(0);
 74: }

 76: PetscErrorCode STSetUp_Sinvert(ST st)
 77: {
 78:   PetscInt       k,nc,nmat=st->nmat;
 79:   PetscScalar    *coeffs=NULL;

 81:   if (nmat>1) STSetWorkVecs(st,1);
 82:   /* if the user did not set the shift, use the target value */
 83:   if (!st->sigma_set) st->sigma = st->defsigma;
 84:   if (nmat>2) {  /* set-up matrices for polynomial eigenproblems */
 85:     if (st->transform) {
 86:       nc = (nmat*(nmat+1))/2;
 87:       PetscMalloc1(nc,&coeffs);
 88:       /* Compute coeffs */
 89:       STCoeffs_Monomial(st,coeffs);
 90:       /* T[0] = A_n */
 91:       k = nmat-1;
 92:       PetscObjectReference((PetscObject)st->A[k]);
 93:       MatDestroy(&st->T[0]);
 94:       st->T[0] = st->A[k];
 95:       for (k=1;k<nmat;k++) 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),PETSC_FALSE,&st->T[k]);
 96:       PetscFree(coeffs);
 97:       PetscObjectReference((PetscObject)st->T[nmat-1]);
 98:       MatDestroy(&st->P);
 99:       st->P = st->T[nmat-1];
100:       if (st->Psplit) {  /* build custom preconditioner from the split matrices */
101:         STMatMAXPY_Private(st,st->sigma,0.0,0,coeffs?coeffs+((nmat-1)*nmat)/2:NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat);
102:       }
103:       ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
104:     } else {
105:       for (k=0;k<nmat;k++) {
106:         PetscObjectReference((PetscObject)st->A[k]);
107:         MatDestroy(&st->T[k]);
108:         st->T[k] = st->A[k];
109:       }
110:     }
111:   }
112:   if (st->P) KSPSetUp(st->ksp);
113:   PetscFunctionReturn(0);
114: }

116: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
117: {
118:   PetscInt       nmat=PetscMax(st->nmat,2),k,nc;
119:   PetscScalar    *coeffs=NULL;

121:   if (st->transform) {
122:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
123:       nc = (nmat*(nmat+1))/2;
124:       PetscMalloc1(nc,&coeffs);
125:       /* Compute coeffs */
126:       STCoeffs_Monomial(st,coeffs);
127:     }
128:     for (k=1;k<nmat;k++) 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,PETSC_FALSE,&st->T[k]);
129:     if (st->matmode == ST_MATMODE_COPY && nmat>2) PetscFree(coeffs);
130:     if (st->P!=st->T[nmat-1]) {
131:       PetscObjectReference((PetscObject)st->T[nmat-1]);
132:       MatDestroy(&st->P);
133:       st->P = st->T[nmat-1];
134:     }
135:     if (st->Psplit) {  /* build custom preconditioner from the split matrices */
136:       STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,0,coeffs?coeffs+((nmat-1)*nmat)/2:NULL,PETSC_FALSE,PETSC_TRUE,&st->Pmat);
137:     }
138:   }
139:   if (st->P) {
140:     ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
141:     KSPSetUp(st->ksp);
142:   }
143:   PetscFunctionReturn(0);
144: }

146: SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
147: {
148:   st->usesksp = PETSC_TRUE;

150:   st->ops->apply           = STApply_Generic;
151:   st->ops->applytrans      = STApplyTranspose_Generic;
152:   st->ops->backtransform   = STBackTransform_Sinvert;
153:   st->ops->setshift        = STSetShift_Sinvert;
154:   st->ops->getbilinearform = STGetBilinearForm_Default;
155:   st->ops->setup           = STSetUp_Sinvert;
156:   st->ops->computeoperator = STComputeOperator_Sinvert;
157:   st->ops->postsolve       = STPostSolve_Sinvert;
158:   st->ops->checknullspace  = STCheckNullSpace_Default;
159:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
160:   PetscFunctionReturn(0);
161: }