Actual source code: sinvert.c

slepc-3.23.0 2025-03-29
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: /*
 17:    Special STApply() for the BSE structured matrix

 19:        H = [ R  C; -C^H -R^T ].

 21:    Assumes that H is a MATNEST and x,y are VECNEST.

 23:    The bottom part of the input vector x2 is computed as
 24:    either conj(x1) or -conj(x1), where the sign is given by
 25:    s in the context SlepcMatStruct.
 26: */
 27: static PetscErrorCode STApply_Sinvert_BSE(ST st,Vec x,Vec y)
 28: {
 29:   Vec            x1,x2;
 30:   PetscContainer container;
 31:   SlepcMatStruct mctx;

 33:   PetscFunctionBegin;
 34:   PetscCall(PetscObjectQuery((PetscObject)st->A[0],"SlepcMatStruct",(PetscObject*)&container));
 35:   PetscCall(PetscContainerGetPointer(container,(void**)&mctx));
 36:   PetscCall(VecNestGetSubVec(x,0,&x1));
 37:   PetscCall(VecNestGetSubVec(x,1,&x2));

 39:   /* x2 = +/-conj(x1) */
 40:   PetscCall(VecCopy(x1,x2));
 41:   PetscCall(VecConjugate(x2));
 42:   if (mctx->s==-1.0) PetscCall(VecScale(x2,-1.0));

 44:   PetscCall(STApply_Generic(st,x,y));
 45:   PetscFunctionReturn(PETSC_SUCCESS);
 46: }

 48: static PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 49: {
 50:   PetscInt    j;
 51: #if !defined(PETSC_USE_COMPLEX)
 52:   PetscScalar t;
 53: #endif

 55:   PetscFunctionBegin;
 56: #if !defined(PETSC_USE_COMPLEX)
 57:   for (j=0;j<n;j++) {
 58:     if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
 59:     else {
 60:       t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
 61:       eigr[j] = eigr[j] / t + st->sigma;
 62:       eigi[j] = - eigi[j] / t;
 63:     }
 64:   }
 65: #else
 66:   for (j=0;j<n;j++) {
 67:     eigr[j] = 1.0 / eigr[j] + st->sigma;
 68:   }
 69: #endif
 70:   PetscFunctionReturn(PETSC_SUCCESS);
 71: }

 73: static PetscErrorCode STPostSolve_Sinvert(ST st)
 74: {
 75:   PetscFunctionBegin;
 76:   if (st->matmode == ST_MATMODE_INPLACE) {
 77:     if (st->nmat>1) PetscCall(MatAXPY(st->A[0],st->sigma,st->A[1],st->str));
 78:     else PetscCall(MatShift(st->A[0],st->sigma));
 79:     st->Astate[0] = ((PetscObject)st->A[0])->state;
 80:     st->state   = ST_STATE_INITIAL;
 81:     st->opready = PETSC_FALSE;
 82:   }
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

 86: /*
 87:    Operator (sinvert):
 88:                Op               P         M
 89:    if nmat=1:  (A-sI)^-1        A-sI      NULL
 90:    if nmat=2:  (A-sB)^-1 B      A-sB      B
 91: */
 92: static PetscErrorCode STComputeOperator_Sinvert(ST st)
 93: {
 94:   PetscFunctionBegin;
 95:   /* if the user did not set the shift, use the target value */
 96:   if (!st->sigma_set) st->sigma = st->defsigma;
 97:   PetscCall(PetscObjectReference((PetscObject)st->A[1]));
 98:   PetscCall(MatDestroy(&st->T[0]));
 99:   st->T[0] = st->A[1];
100:   PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]));
101:   PetscCall(PetscObjectReference((PetscObject)st->T[1]));
102:   PetscCall(MatDestroy(&st->P));
103:   st->P = st->T[1];
104:   st->M = (st->nmat>1)? st->T[0]: NULL;
105:   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
106:     PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
107:   }
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: static PetscErrorCode STSetUp_Sinvert(ST st)
112: {
113:   PetscInt       k,nc,nmat=st->nmat,nsub;
114:   PetscScalar    *coeffs=NULL;
115:   PetscBool      islu;
116:   KSP            *subksp;
117:   PC             pc,pcsub;
118:   Mat            A00;
119:   MatType        type;
120:   PetscBool      flg;
121:   char           str[64];

123:   PetscFunctionBegin;
124:   if (nmat>1) PetscCall(STSetWorkVecs(st,1));
125:   /* if the user did not set the shift, use the target value */
126:   if (!st->sigma_set) st->sigma = st->defsigma;
127:   if (nmat>2) {  /* set-up matrices for polynomial eigenproblems */
128:     if (st->transform) {
129:       nc = (nmat*(nmat+1))/2;
130:       PetscCall(PetscMalloc1(nc,&coeffs));
131:       /* Compute coeffs */
132:       PetscCall(STCoeffs_Monomial(st,coeffs));
133:       /* T[0] = A_n */
134:       k = nmat-1;
135:       PetscCall(PetscObjectReference((PetscObject)st->A[k]));
136:       PetscCall(MatDestroy(&st->T[0]));
137:       st->T[0] = st->A[k];
138:       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]));
139:       PetscCall(PetscFree(coeffs));
140:       PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
141:       PetscCall(MatDestroy(&st->P));
142:       st->P = st->T[nmat-1];
143:       if (st->Psplit) {  /* build custom preconditioner from the split matrices */
144:         PetscCall(STMatMAXPY_Private(st,st->sigma,0.0,0,coeffs?coeffs+((nmat-1)*nmat)/2:NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
145:       }
146:       PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
147:     } else {
148:       for (k=0;k<nmat;k++) {
149:         PetscCall(PetscObjectReference((PetscObject)st->A[k]));
150:         PetscCall(MatDestroy(&st->T[k]));
151:         st->T[k] = st->A[k];
152:       }
153:     }
154:   }
155:   if (st->structured) {
156:     /* ./ex55 -st_type sinvert -eps_target 0 -eps_target_magnitude
157:               -st_ksp_type preonly -st_pc_type fieldsplit
158:               -st_fieldsplit_0_ksp_type preonly -st_fieldsplit_0_pc_type lu
159:               -st_fieldsplit_1_ksp_type preonly -st_fieldsplit_1_pc_type lu
160:               -st_pc_fieldsplit_type schur -st_pc_fieldsplit_schur_fact_type full
161:               -st_pc_fieldsplit_schur_precondition full */
162:     PetscCall(KSPGetPC(st->ksp,&pc));
163:     PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCLU,&islu));
164:     if (islu) {  /* assume PCLU means the user has not set any options */
165:       PetscCall(KSPSetType(st->ksp,KSPPREONLY));
166:       PetscCall(PCSetType(pc,PCFIELDSPLIT));
167:       PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR));
168:       PetscCall(PCFieldSplitSetSchurFactType(pc,PC_FIELDSPLIT_SCHUR_FACT_FULL));
169:       PetscCall(PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_FULL,NULL));
170:       /* hack to set Mat type of Schur complement equal to A00, assumes default prefixes */
171:       PetscCall(PetscOptionsHasName(((PetscObject)st)->options,((PetscObject)st)->prefix,"-st_fieldsplit_1_explicit_operator_mat_type",&flg));
172:       if (!flg) {
173:         PetscCall(MatNestGetSubMat(st->A[0],0,0,&A00));
174:         PetscCall(MatGetType(A00,&type));
175:         PetscCall(PetscSNPrintf(str,sizeof(str),"-%sst_fieldsplit_1_explicit_operator_mat_type %s",((PetscObject)st)->prefix?((PetscObject)st)->prefix:"",type));
176:         PetscCall(PetscOptionsInsertString(((PetscObject)st)->options,str));
177:       }
178:       /* set preonly+lu on block solvers */
179:       PetscCall(KSPSetUp(st->ksp));
180:       PetscCall(PCFieldSplitGetSubKSP(pc,&nsub,&subksp));
181:       PetscCall(KSPSetType(subksp[0],KSPPREONLY));
182:       PetscCall(KSPGetPC(subksp[0],&pcsub));
183:       PetscCall(PCSetType(pcsub,PCLU));
184:       PetscCall(KSPSetType(subksp[1],KSPPREONLY));
185:       PetscCall(KSPGetPC(subksp[1],&pcsub));
186:       PetscCall(PCSetType(pcsub,PCLU));
187:       PetscCall(PetscFree(subksp));
188:     }
189:     st->ops->apply = STApply_Sinvert_BSE;
190:   } else {
191:     if (st->P) PetscCall(KSPSetUp(st->ksp));
192:   }
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: static PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
197: {
198:   PetscInt       nmat=PetscMax(st->nmat,2),k,nc;
199:   PetscScalar    *coeffs=NULL;

201:   PetscFunctionBegin;
202:   if (st->transform) {
203:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
204:       nc = (nmat*(nmat+1))/2;
205:       PetscCall(PetscMalloc1(nc,&coeffs));
206:       /* Compute coeffs */
207:       PetscCall(STCoeffs_Monomial(st,coeffs));
208:     }
209:     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]));
210:     if (st->matmode == ST_MATMODE_COPY && nmat>2) PetscCall(PetscFree(coeffs));
211:     if (st->P!=st->T[nmat-1]) {
212:       PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
213:       PetscCall(MatDestroy(&st->P));
214:       st->P = st->T[nmat-1];
215:     }
216:     if (st->Psplit) {  /* build custom preconditioner from the split matrices */
217:       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));
218:     }
219:   }
220:   if (st->P) {
221:     PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
222:     PetscCall(KSPSetUp(st->ksp));
223:   }
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
228: {
229:   PetscFunctionBegin;
230:   st->usesksp = PETSC_TRUE;

232:   st->ops->apply           = STApply_Generic;
233:   st->ops->applytrans      = STApplyTranspose_Generic;
234:   st->ops->applyhermtrans  = STApplyHermitianTranspose_Generic;
235:   st->ops->backtransform   = STBackTransform_Sinvert;
236:   st->ops->setshift        = STSetShift_Sinvert;
237:   st->ops->getbilinearform = STGetBilinearForm_Default;
238:   st->ops->setup           = STSetUp_Sinvert;
239:   st->ops->computeoperator = STComputeOperator_Sinvert;
240:   st->ops->postsolve       = STPostSolve_Sinvert;
241:   st->ops->checknullspace  = STCheckNullSpace_Default;
242:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }