Actual source code: stsolve.c

slepc-3.21.0 2024-03-30
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:    ST interface routines, callable by users
 12: */

 14: #include <slepc/private/stimpl.h>

 16: PetscErrorCode STApply_Generic(ST st,Vec x,Vec y)
 17: {
 18:   PetscFunctionBegin;
 19:   if (st->M && st->P) {
 20:     PetscCall(MatMult(st->M,x,st->work[0]));
 21:     PetscCall(STMatSolve(st,st->work[0],y));
 22:   } else if (st->M) PetscCall(MatMult(st->M,x,y));
 23:   else PetscCall(STMatSolve(st,x,y));
 24:   PetscFunctionReturn(PETSC_SUCCESS);
 25: }

 27: /*@
 28:    STApply - Applies the spectral transformation operator to a vector, for
 29:    instance (A - sB)^-1 B in the case of the shift-and-invert transformation
 30:    and generalized eigenproblem.

 32:    Collective

 34:    Input Parameters:
 35: +  st - the spectral transformation context
 36: -  x  - input vector

 38:    Output Parameter:
 39: .  y - output vector

 41:    Level: developer

 43: .seealso: STApplyTranspose(), STApplyHermitianTranspose()
 44: @*/
 45: PetscErrorCode STApply(ST st,Vec x,Vec y)
 46: {
 47:   Mat            Op;

 49:   PetscFunctionBegin;
 54:   STCheckMatrices(st,1);
 55:   PetscCheck(x!=y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
 56:   PetscCall(VecSetErrorIfLocked(y,3));
 57:   PetscCall(STGetOperator_Private(st,&Op));
 58:   PetscCall(MatMult(Op,x,y));
 59:   PetscFunctionReturn(PETSC_SUCCESS);
 60: }

 62: PetscErrorCode STApplyMat_Generic(ST st,Mat B,Mat C)
 63: {
 64:   Mat            work;

 66:   PetscFunctionBegin;
 67:   if (st->M && st->P) {
 68:     PetscCall(MatMatMult(st->M,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&work));
 69:     PetscCall(STMatMatSolve(st,work,C));
 70:     PetscCall(MatDestroy(&work));
 71:   } else if (st->M) PetscCall(MatMatMult(st->M,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C));
 72:   else PetscCall(STMatMatSolve(st,B,C));
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: /*@
 77:    STApplyMat - Applies the spectral transformation operator to a matrix, for
 78:    instance (A - sB)^-1 B in the case of the shift-and-invert transformation
 79:    and generalized eigenproblem.

 81:    Collective

 83:    Input Parameters:
 84: +  st - the spectral transformation context
 85: -  X  - input matrix

 87:    Output Parameter:
 88: .  Y - output matrix

 90:    Level: developer

 92: .seealso: STApply()
 93: @*/
 94: PetscErrorCode STApplyMat(ST st,Mat X,Mat Y)
 95: {
 96:   PetscFunctionBegin;
101:   STCheckMatrices(st,1);
102:   PetscCheck(X!=Y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"X and Y must be different matrices");
103:   PetscUseTypeMethod(st,applymat,X,Y);
104:   PetscFunctionReturn(PETSC_SUCCESS);
105: }

107: PetscErrorCode STApplyTranspose_Generic(ST st,Vec x,Vec y)
108: {
109:   PetscFunctionBegin;
110:   if (st->M && st->P) {
111:     PetscCall(STMatSolveTranspose(st,x,st->work[0]));
112:     PetscCall(MatMultTranspose(st->M,st->work[0],y));
113:   } else if (st->M) PetscCall(MatMultTranspose(st->M,x,y));
114:   else PetscCall(STMatSolveTranspose(st,x,y));
115:   PetscFunctionReturn(PETSC_SUCCESS);
116: }

118: /*@
119:    STApplyTranspose - Applies the transpose of the operator to a vector, for
120:    instance B^T(A - sB)^-T in the case of the shift-and-invert transformation
121:    and generalized eigenproblem.

123:    Collective

125:    Input Parameters:
126: +  st - the spectral transformation context
127: -  x  - input vector

129:    Output Parameter:
130: .  y - output vector

132:    Level: developer

134: .seealso: STApply(), STApplyHermitianTranspose()
135: @*/
136: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
137: {
138:   Mat            Op;

140:   PetscFunctionBegin;
145:   STCheckMatrices(st,1);
146:   PetscCheck(x!=y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
147:   PetscCall(VecSetErrorIfLocked(y,3));
148:   PetscCall(STGetOperator_Private(st,&Op));
149:   PetscCall(MatMultTranspose(Op,x,y));
150:   PetscFunctionReturn(PETSC_SUCCESS);
151: }

153: PetscErrorCode STApplyHermitianTranspose_Generic(ST st,Vec x,Vec y)
154: {
155:   PetscFunctionBegin;
156:   if (st->M && st->P) {
157:     PetscCall(STMatSolveHermitianTranspose(st,x,st->work[0]));
158:     PetscCall(MatMultHermitianTranspose(st->M,st->work[0],y));
159:   } else if (st->M) PetscCall(MatMultHermitianTranspose(st->M,x,y));
160:   else PetscCall(STMatSolveHermitianTranspose(st,x,y));
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }

164: /*@
165:    STApplyHermitianTranspose - Applies the hermitian-transpose of the operator
166:    to a vector, for instance B^H(A - sB)^-H in the case of the shift-and-invert
167:    transformation and generalized eigenproblem.

169:    Collective

171:    Input Parameters:
172: +  st - the spectral transformation context
173: -  x  - input vector

175:    Output Parameter:
176: .  y - output vector

178:    Note:
179:    Currently implemented via STApplyTranspose() with appropriate conjugation.

181:    Level: developer

183: .seealso: STApply(), STApplyTranspose()
184: @*/
185: PetscErrorCode STApplyHermitianTranspose(ST st,Vec x,Vec y)
186: {
187:   Mat            Op;

189:   PetscFunctionBegin;
194:   STCheckMatrices(st,1);
195:   PetscCheck(x!=y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
196:   PetscCall(VecSetErrorIfLocked(y,3));
197:   PetscCall(STGetOperator_Private(st,&Op));
198:   PetscCall(MatMultHermitianTranspose(Op,x,y));
199:   PetscFunctionReturn(PETSC_SUCCESS);
200: }

202: /*@
203:    STGetBilinearForm - Returns the matrix used in the bilinear form with a
204:    generalized problem with semi-definite B.

206:    Logically Collective

208:    Input Parameters:
209: .  st - the spectral transformation context

211:    Output Parameter:
212: .  B - output matrix

214:    Notes:
215:    The output matrix B must be destroyed after use. It will be NULL in
216:    case of standard eigenproblems.

218:    Level: developer

220: .seealso: BVSetMatrix()
221: @*/
222: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
223: {
224:   PetscFunctionBegin;
227:   PetscAssertPointer(B,2);
228:   STCheckMatrices(st,1);
229:   PetscUseTypeMethod(st,getbilinearform,B);
230:   PetscFunctionReturn(PETSC_SUCCESS);
231: }

233: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
234: {
235:   PetscFunctionBegin;
236:   if (st->nmat==1) *B = NULL;
237:   else {
238:     *B = st->A[1];
239:     PetscCall(PetscObjectReference((PetscObject)*B));
240:   }
241:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

244: static PetscErrorCode MatMult_STOperator(Mat Op,Vec x,Vec y)
245: {
246:   ST             st;

248:   PetscFunctionBegin;
249:   PetscCall(MatShellGetContext(Op,&st));
250:   PetscCall(STSetUp(st));
251:   PetscCall(PetscLogEventBegin(ST_Apply,st,x,y,0));
252:   if (st->D) { /* with balancing */
253:     PetscCall(VecPointwiseDivide(st->wb,x,st->D));
254:     PetscUseTypeMethod(st,apply,st->wb,y);
255:     PetscCall(VecPointwiseMult(y,y,st->D));
256:   } else PetscUseTypeMethod(st,apply,x,y);
257:   PetscCall(PetscLogEventEnd(ST_Apply,st,x,y,0));
258:   PetscFunctionReturn(PETSC_SUCCESS);
259: }

261: static PetscErrorCode MatMultTranspose_STOperator(Mat Op,Vec x,Vec y)
262: {
263:   ST             st;

265:   PetscFunctionBegin;
266:   PetscCall(MatShellGetContext(Op,&st));
267:   PetscCall(STSetUp(st));
268:   PetscCall(PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0));
269:   if (st->D) { /* with balancing */
270:     PetscCall(VecPointwiseMult(st->wb,x,st->D));
271:     PetscUseTypeMethod(st,applytrans,st->wb,y);
272:     PetscCall(VecPointwiseDivide(y,y,st->D));
273:   } else PetscUseTypeMethod(st,applytrans,x,y);
274:   PetscCall(PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0));
275:   PetscFunctionReturn(PETSC_SUCCESS);
276: }

278: #if defined(PETSC_USE_COMPLEX)
279: static PetscErrorCode MatMultHermitianTranspose_STOperator(Mat Op,Vec x,Vec y)
280: {
281:   ST             st;

283:   PetscFunctionBegin;
284:   PetscCall(MatShellGetContext(Op,&st));
285:   PetscCall(STSetUp(st));
286:   PetscCall(PetscLogEventBegin(ST_ApplyHermitianTranspose,st,x,y,0));
287:   if (st->ops->applyhermtrans) {
288:     if (st->D) { /* with balancing */
289:       PetscCall(VecConjugate(st->D));
290:       PetscCall(VecPointwiseMult(st->wb,x,st->D));
291:       PetscUseTypeMethod(st,applyhermtrans,st->wb,y);
292:       PetscCall(VecPointwiseDivide(y,y,st->D));
293:       PetscCall(VecConjugate(st->D));
294:     } else PetscUseTypeMethod(st,applyhermtrans,x,y);
295:   } else {
296:     if (!st->wht) PetscCall(MatCreateVecs(st->A[0],&st->wht,NULL));
297:     PetscCall(VecCopy(x,st->wht));
298:     PetscCall(VecConjugate(st->wht));
299:     if (st->D) { /* with balancing */
300:       PetscCall(VecPointwiseMult(st->wb,st->wht,st->D));
301:       PetscUseTypeMethod(st,applytrans,st->wb,y);
302:       PetscCall(VecPointwiseDivide(y,y,st->D));
303:     } else PetscUseTypeMethod(st,applytrans,st->wht,y);
304:     PetscCall(VecConjugate(y));
305:   }
306:   PetscCall(PetscLogEventEnd(ST_ApplyHermitianTranspose,st,x,y,0));
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }
309: #endif

311: static PetscErrorCode MatMatMult_STOperator(Mat Op,Mat B,Mat C,void *ctx)
312: {
313:   ST             st;

315:   PetscFunctionBegin;
316:   PetscCall(MatShellGetContext(Op,&st));
317:   PetscCall(STSetUp(st));
318:   PetscCall(PetscLogEventBegin(ST_Apply,st,B,C,0));
319:   PetscCall(STApplyMat_Generic(st,B,C));
320:   PetscCall(PetscLogEventEnd(ST_Apply,st,B,C,0));
321:   PetscFunctionReturn(PETSC_SUCCESS);
322: }

324: PetscErrorCode STGetOperator_Private(ST st,Mat *Op)
325: {
326:   PetscInt       m,n,M,N;
327:   Vec            v;
328:   VecType        vtype;

330:   PetscFunctionBegin;
331:   if (!st->Op) {
332:     if (Op) *Op = NULL;
333:     /* create the shell matrix */
334:     PetscCall(MatGetLocalSize(st->A[0],&m,&n));
335:     PetscCall(MatGetSize(st->A[0],&M,&N));
336:     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)st),m,n,M,N,st,&st->Op));
337:     PetscCall(MatShellSetOperation(st->Op,MATOP_MULT,(void(*)(void))MatMult_STOperator));
338:     PetscCall(MatShellSetOperation(st->Op,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_STOperator));
339: #if defined(PETSC_USE_COMPLEX)
340:     PetscCall(MatShellSetOperation(st->Op,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultHermitianTranspose_STOperator));
341: #else
342:     PetscCall(MatShellSetOperation(st->Op,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultTranspose_STOperator));
343: #endif
344:     if (!st->D && st->ops->apply==STApply_Generic) {
345:       PetscCall(MatShellSetMatProductOperation(st->Op,MATPRODUCT_AB,NULL,MatMatMult_STOperator,NULL,MATDENSE,MATDENSE));
346:       PetscCall(MatShellSetMatProductOperation(st->Op,MATPRODUCT_AB,NULL,MatMatMult_STOperator,NULL,MATDENSECUDA,MATDENSECUDA));
347:     }
348:     /* make sure the shell matrix generates a vector of the same type as the problem matrices */
349:     PetscCall(MatCreateVecs(st->A[0],&v,NULL));
350:     PetscCall(VecGetType(v,&vtype));
351:     PetscCall(MatShellSetVecType(st->Op,vtype));
352:     PetscCall(VecDestroy(&v));
353:     /* build the operator matrices */
354:     PetscCall(STComputeOperator(st));
355:   }
356:   if (Op) *Op = st->Op;
357:   PetscFunctionReturn(PETSC_SUCCESS);
358: }

360: /*@
361:    STGetOperator - Returns a shell matrix that represents the operator of the
362:    spectral transformation.

364:    Collective

366:    Input Parameter:
367: .  st - the spectral transformation context

369:    Output Parameter:
370: .  Op - operator matrix

372:    Notes:
373:    The operator is defined in linear eigenproblems only, not in polynomial ones,
374:    so the call will fail if more than 2 matrices were passed in STSetMatrices().

376:    The returned shell matrix is essentially a wrapper to the STApply() and
377:    STApplyTranspose() operations. The operator can often be expressed as

379: $     Op = D*inv(K)*M*inv(D)

381:    where D is the balancing matrix, and M and K are two matrices corresponding
382:    to the numerator and denominator for spectral transformations that represent
383:    a rational matrix function. In the case of STSHELL, the inner part inv(K)*M
384:    is replaced by the user-provided operation from STShellSetApply().

386:    The preconditioner matrix K typically depends on the value of the shift, and
387:    its inverse is handled via an internal KSP object. Normal usage does not
388:    require explicitly calling STGetOperator(), but it can be used to force the
389:    creation of K and M, and then K is passed to the KSP. This is useful for
390:    setting options associated with the PCFactor (to set MUMPS options, for instance).

392:    The returned matrix must NOT be destroyed by the user. Instead, when no
393:    longer needed it must be returned with STRestoreOperator(). In particular,
394:    this is required before modifying the ST matrices or the shift.

396:    A NULL pointer can be passed in Op in case the matrix is not required but we
397:    want to force its creation. In this case, STRestoreOperator() should not be
398:    called.

400:    Level: advanced

402: .seealso: STApply(), STApplyTranspose(), STSetBalanceMatrix(), STShellSetApply(),
403:           STGetKSP(), STSetShift(), STRestoreOperator(), STSetMatrices()
404: @*/
405: PetscErrorCode STGetOperator(ST st,Mat *Op)
406: {
407:   PetscFunctionBegin;
410:   STCheckMatrices(st,1);
411:   STCheckNotSeized(st,1);
412:   PetscCheck(st->nmat<=2,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"The operator is not defined in polynomial eigenproblems");
413:   PetscCall(STGetOperator_Private(st,Op));
414:   if (Op) st->opseized = PETSC_TRUE;
415:   PetscFunctionReturn(PETSC_SUCCESS);
416: }

418: /*@
419:    STRestoreOperator - Restore the previously seized operator matrix.

421:    Logically Collective

423:    Input Parameters:
424: +  st - the spectral transformation context
425: -  Op - operator matrix

427:    Notes:
428:    The arguments must match the corresponding call to STGetOperator().

430:    Level: advanced

432: .seealso: STGetOperator()
433: @*/
434: PetscErrorCode STRestoreOperator(ST st,Mat *Op)
435: {
436:   PetscFunctionBegin;
438:   PetscAssertPointer(Op,2);
440:   PetscCheck(st->opseized,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Must be called after STGetOperator()");
441:   *Op = NULL;
442:   st->opseized = PETSC_FALSE;
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: /*
447:    STComputeOperator - Computes the matrices that constitute the operator

449:       Op = D*inv(K)*M*inv(D).

451:    K and M are computed here (D is user-provided) from the system matrices
452:    and the shift sigma (whenever these are changed, this function recomputes
453:    K and M). This is used only in linear eigenproblems (nmat<3).

455:    K is the "preconditioner matrix": it is the denominator in rational operators,
456:    e.g. (A-sigma*B) in shift-and-invert. In non-rational transformations such
457:    as STFILTER, K=NULL which means identity. After computing K, it is passed to
458:    the internal KSP object via KSPSetOperators.

460:    M is the numerator in rational operators. If unused it is set to NULL (e.g.
461:    in STPRECOND).

463:    STSHELL does not compute anything here, but sets the flag as if it was ready.
464: */
465: PetscErrorCode STComputeOperator(ST st)
466: {
467:   PC             pc;

469:   PetscFunctionBegin;
472:   if (!st->opready && st->ops->computeoperator) {
473:     PetscCall(PetscInfo(st,"Building the operator matrices\n"));
474:     STCheckMatrices(st,1);
475:     if (!st->T) PetscCall(PetscCalloc1(PetscMax(2,st->nmat),&st->T));
476:     PetscCall(PetscLogEventBegin(ST_ComputeOperator,st,0,0,0));
477:     PetscUseTypeMethod(st,computeoperator);
478:     PetscCall(PetscLogEventEnd(ST_ComputeOperator,st,0,0,0));
479:     if (st->usesksp) {
480:       if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp));
481:       if (st->P) {
482:         PetscCall(STSetDefaultKSP(st));
483:         PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
484:       } else {
485:         /* STPRECOND defaults to PCNONE if st->P is empty */
486:         PetscCall(KSPGetPC(st->ksp,&pc));
487:         PetscCall(PCSetType(pc,PCNONE));
488:       }
489:     }
490:   }
491:   st->opready = PETSC_TRUE;
492:   PetscFunctionReturn(PETSC_SUCCESS);
493: }

495: /*@
496:    STSetUp - Prepares for the use of a spectral transformation.

498:    Collective

500:    Input Parameter:
501: .  st - the spectral transformation context

503:    Level: advanced

505: .seealso: STCreate(), STApply(), STDestroy()
506: @*/
507: PetscErrorCode STSetUp(ST st)
508: {
509:   PetscInt       i,n,k;

511:   PetscFunctionBegin;
514:   STCheckMatrices(st,1);
515:   switch (st->state) {
516:     case ST_STATE_INITIAL:
517:       PetscCall(PetscInfo(st,"Setting up new ST\n"));
518:       if (!((PetscObject)st)->type_name) PetscCall(STSetType(st,STSHIFT));
519:       break;
520:     case ST_STATE_SETUP:
521:       PetscFunctionReturn(PETSC_SUCCESS);
522:     case ST_STATE_UPDATED:
523:       PetscCall(PetscInfo(st,"Setting up updated ST\n"));
524:       break;
525:   }
526:   PetscCall(PetscLogEventBegin(ST_SetUp,st,0,0,0));
527:   if (st->state!=ST_STATE_UPDATED) {
528:     if (!(st->nmat<3 && st->opready)) {
529:       if (st->T) {
530:         for (i=0;i<PetscMax(2,st->nmat);i++) PetscCall(MatDestroy(&st->T[i]));
531:       }
532:       PetscCall(MatDestroy(&st->P));
533:     }
534:   }
535:   if (st->D) {
536:     PetscCall(MatGetLocalSize(st->A[0],NULL,&n));
537:     PetscCall(VecGetLocalSize(st->D,&k));
538:     PetscCheck(n==k,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_SIZ,"Balance matrix has wrong dimension %" PetscInt_FMT " (should be %" PetscInt_FMT ")",k,n);
539:     if (!st->wb) PetscCall(VecDuplicate(st->D,&st->wb));
540:   }
541:   if (st->nmat<3 && st->transform) PetscCall(STComputeOperator(st));
542:   else {
543:     if (!st->T) PetscCall(PetscCalloc1(PetscMax(2,st->nmat),&st->T));
544:   }
545:   PetscTryTypeMethod(st,setup);
546:   st->state = ST_STATE_SETUP;
547:   PetscCall(PetscLogEventEnd(ST_SetUp,st,0,0,0));
548:   PetscFunctionReturn(PETSC_SUCCESS);
549: }

551: /*
552:    Computes coefficients for the transformed polynomial,
553:    and stores the result in argument S.

555:    alpha - value of the parameter of the transformed polynomial
556:    beta - value of the previous shift (only used in inplace mode)
557:    k - index of first matrix included in the computation
558:    coeffs - coefficients of the expansion
559:    initial - true if this is the first time
560:    precond - whether the preconditioner matrix must be computed
561: */
562: PetscErrorCode STMatMAXPY_Private(ST st,PetscScalar alpha,PetscScalar beta,PetscInt k,PetscScalar *coeffs,PetscBool initial,PetscBool precond,Mat *S)
563: {
564:   PetscInt       *matIdx=NULL,nmat,i,ini=-1;
565:   PetscScalar    t=1.0,ta,gamma;
566:   PetscBool      nz=PETSC_FALSE;
567:   Mat            *A=precond?st->Psplit:st->A;
568:   MatStructure   str=precond?st->strp:st->str;

570:   PetscFunctionBegin;
571:   nmat = st->nmat-k;
572:   switch (st->matmode) {
573:   case ST_MATMODE_INPLACE:
574:     PetscCheck(st->nmat<=2,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST_MATMODE_INPLACE not supported for polynomial eigenproblems");
575:     PetscCheck(!precond,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST_MATMODE_INPLACE not supported for split preconditioner");
576:     if (initial) {
577:       PetscCall(PetscObjectReference((PetscObject)A[0]));
578:       *S = A[0];
579:       gamma = alpha;
580:     } else gamma = alpha-beta;
581:     if (gamma != 0.0) {
582:       if (st->nmat>1) PetscCall(MatAXPY(*S,gamma,A[1],str));
583:       else PetscCall(MatShift(*S,gamma));
584:     }
585:     break;
586:   case ST_MATMODE_SHELL:
587:     PetscCheck(!precond,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST_MATMODE_SHELL not supported for split preconditioner");
588:     if (initial) {
589:       if (st->nmat>2) {
590:         PetscCall(PetscMalloc1(nmat,&matIdx));
591:         for (i=0;i<nmat;i++) matIdx[i] = k+i;
592:       }
593:       PetscCall(STMatShellCreate(st,alpha,nmat,matIdx,coeffs,S));
594:       if (st->nmat>2) PetscCall(PetscFree(matIdx));
595:     } else PetscCall(STMatShellShift(*S,alpha));
596:     break;
597:   case ST_MATMODE_COPY:
598:     if (coeffs) {
599:       for (i=0;i<nmat && ini==-1;i++) {
600:         if (coeffs[i]!=0.0) ini = i;
601:         else t *= alpha;
602:       }
603:       if (coeffs[ini] != 1.0) nz = PETSC_TRUE;
604:       for (i=ini+1;i<nmat&&!nz;i++) if (coeffs[i]!=0.0) nz = PETSC_TRUE;
605:     } else { nz = PETSC_TRUE; ini = 0; }
606:     if ((alpha == 0.0 || !nz) && t==1.0) {
607:       PetscCall(PetscObjectReference((PetscObject)A[k+ini]));
608:       PetscCall(MatDestroy(S));
609:       *S = A[k+ini];
610:     } else {
611:       if (*S && *S!=A[k+ini]) {
612:         PetscCall(MatSetOption(*S,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
613:         PetscCall(MatCopy(A[k+ini],*S,DIFFERENT_NONZERO_PATTERN));
614:       } else {
615:         PetscCall(MatDestroy(S));
616:         PetscCall(MatDuplicate(A[k+ini],MAT_COPY_VALUES,S));
617:         PetscCall(MatSetOption(*S,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
618:       }
619:       if (coeffs && coeffs[ini]!=1.0) PetscCall(MatScale(*S,coeffs[ini]));
620:       for (i=ini+k+1;i<PetscMax(2,st->nmat);i++) {
621:         t *= alpha;
622:         ta = t;
623:         if (coeffs) ta *= coeffs[i-k];
624:         if (ta!=0.0) {
625:           if (st->nmat>1) PetscCall(MatAXPY(*S,ta,A[i],str));
626:           else PetscCall(MatShift(*S,ta));
627:         }
628:       }
629:     }
630:   }
631:   PetscCall(MatSetOption(*S,MAT_SYMMETRIC,st->asymm));
632:   PetscCall(MatSetOption(*S,MAT_HERMITIAN,(PetscImaginaryPart(st->sigma)==0.0)?st->aherm:PETSC_FALSE));
633:   PetscFunctionReturn(PETSC_SUCCESS);
634: }

636: /*
637:    Computes the values of the coefficients required by STMatMAXPY_Private
638:    for the case of monomial basis.
639: */
640: PetscErrorCode STCoeffs_Monomial(ST st, PetscScalar *coeffs)
641: {
642:   PetscInt  k,i,ini,inip;

644:   PetscFunctionBegin;
645:   /* Compute binomial coefficients */
646:   ini = (st->nmat*(st->nmat-1))/2;
647:   for (i=0;i<st->nmat;i++) coeffs[ini+i]=1.0;
648:   for (k=st->nmat-1;k>=1;k--) {
649:     inip = ini+1;
650:     ini = (k*(k-1))/2;
651:     coeffs[ini] = 1.0;
652:     for (i=1;i<k;i++) coeffs[ini+i] = coeffs[ini+i-1]+coeffs[inip+i-1];
653:   }
654:   PetscFunctionReturn(PETSC_SUCCESS);
655: }

657: /*@
658:    STPostSolve - Optional post-solve phase, intended for any actions that must
659:    be performed on the ST object after the eigensolver has finished.

661:    Collective

663:    Input Parameters:
664: .  st  - the spectral transformation context

666:    Level: developer

668: .seealso: EPSSolve()
669: @*/
670: PetscErrorCode STPostSolve(ST st)
671: {
672:   PetscFunctionBegin;
675:   PetscTryTypeMethod(st,postsolve);
676:   PetscFunctionReturn(PETSC_SUCCESS);
677: }

679: /*@
680:    STBackTransform - Back-transformation phase, intended for
681:    spectral transformations which require to transform the computed
682:    eigenvalues back to the original eigenvalue problem.

684:    Not Collective

686:    Input Parameters:
687: +  st   - the spectral transformation context
688: .  n    - number of eigenvalues
689: .  eigr - real part of a computed eigenvalues
690: -  eigi - imaginary part of a computed eigenvalues

692:    Level: developer

694: .seealso: STIsInjective()
695: @*/
696: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
697: {
698:   PetscFunctionBegin;
701:   PetscTryTypeMethod(st,backtransform,n,eigr,eigi);
702:   PetscFunctionReturn(PETSC_SUCCESS);
703: }

705: /*@
706:    STIsInjective - Ask if this spectral transformation is injective or not
707:    (that is, if it corresponds to a one-to-one mapping). If not, then it
708:    does not make sense to call STBackTransform().

710:    Not Collective

712:    Input Parameter:
713: .  st   - the spectral transformation context

715:    Output Parameter:
716: .  is - the answer

718:    Level: developer

720: .seealso: STBackTransform()
721: @*/
722: PetscErrorCode STIsInjective(ST st,PetscBool* is)
723: {
724:   PetscBool      shell;

726:   PetscFunctionBegin;
729:   PetscAssertPointer(is,2);

731:   PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&shell));
732:   if (shell) PetscCall(STIsInjective_Shell(st,is));
733:   else *is = st->ops->backtransform? PETSC_TRUE: PETSC_FALSE;
734:   PetscFunctionReturn(PETSC_SUCCESS);
735: }

737: /*@
738:    STMatSetUp - Build the preconditioner matrix used in STMatSolve().

740:    Collective

742:    Input Parameters:
743: +  st     - the spectral transformation context
744: .  sigma  - the shift
745: -  coeffs - the coefficients (may be NULL)

747:    Note:
748:    This function is not intended to be called by end users, but by SLEPc
749:    solvers that use ST. It builds matrix st->P as follows, then calls KSPSetUp().
750: .vb
751:     If (coeffs)  st->P = Sum_{i=0..nmat-1} coeffs[i]*sigma^i*A_i
752:     else         st->P = Sum_{i=0..nmat-1} sigma^i*A_i
753: .ve

755:    Level: developer

757: .seealso: STMatSolve()
758: @*/
759: PetscErrorCode STMatSetUp(ST st,PetscScalar sigma,PetscScalar *coeffs)
760: {
761:   PetscFunctionBegin;
764:   STCheckMatrices(st,1);

766:   PetscCall(PetscLogEventBegin(ST_MatSetUp,st,0,0,0));
767:   PetscCall(STMatMAXPY_Private(st,sigma,0.0,0,coeffs,PETSC_TRUE,PETSC_FALSE,&st->P));
768:   if (st->Psplit) PetscCall(STMatMAXPY_Private(st,sigma,0.0,0,coeffs,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
769:   PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
770:   PetscCall(KSPSetUp(st->ksp));
771:   PetscCall(PetscLogEventEnd(ST_MatSetUp,st,0,0,0));
772:   PetscFunctionReturn(PETSC_SUCCESS);
773: }

775: /*@
776:    STSetWorkVecs - Sets a number of work vectors into the ST object.

778:    Collective

780:    Input Parameters:
781: +  st - the spectral transformation context
782: -  nw - number of work vectors to allocate

784:    Developer Notes:
785:    This is SLEPC_EXTERN because it may be required by shell STs.

787:    Level: developer

789: .seealso: STMatCreateVecs()
790: @*/
791: PetscErrorCode STSetWorkVecs(ST st,PetscInt nw)
792: {
793:   PetscInt       i;

795:   PetscFunctionBegin;
798:   PetscCheck(nw>0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"nw must be > 0: nw = %" PetscInt_FMT,nw);
799:   if (st->nwork < nw) {
800:     PetscCall(VecDestroyVecs(st->nwork,&st->work));
801:     st->nwork = nw;
802:     PetscCall(PetscMalloc1(nw,&st->work));
803:     for (i=0;i<nw;i++) PetscCall(STMatCreateVecs(st,&st->work[i],NULL));
804:   }
805:   PetscFunctionReturn(PETSC_SUCCESS);
806: }