Actual source code: sinvert.c

slepc-3.21.1 2024-04-26
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: static 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:   PetscFunctionBegin;
 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:   PetscFunctionReturn(PETSC_SUCCESS);
 39: }

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

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

 79: static PetscErrorCode STSetUp_Sinvert(ST st)
 80: {
 81:   PetscInt       k,nc,nmat=st->nmat;
 82:   PetscScalar    *coeffs=NULL;

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

120: static PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
121: {
122:   PetscInt       nmat=PetscMax(st->nmat,2),k,nc;
123:   PetscScalar    *coeffs=NULL;

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

151: SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
152: {
153:   PetscFunctionBegin;
154:   st->usesksp = PETSC_TRUE;

156:   st->ops->apply           = STApply_Generic;
157:   st->ops->applytrans      = STApplyTranspose_Generic;
158:   st->ops->applyhermtrans  = STApplyHermitianTranspose_Generic;
159:   st->ops->backtransform   = STBackTransform_Sinvert;
160:   st->ops->setshift        = STSetShift_Sinvert;
161:   st->ops->getbilinearform = STGetBilinearForm_Default;
162:   st->ops->setup           = STSetUp_Sinvert;
163:   st->ops->computeoperator = STComputeOperator_Sinvert;
164:   st->ops->postsolve       = STPostSolve_Sinvert;
165:   st->ops->checknullspace  = STCheckNullSpace_Default;
166:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
167:   PetscFunctionReturn(PETSC_SUCCESS);
168: }