Actual source code: stsolve.c

slepc-main 2024-11-15
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:       PetscCall(MatShellSetMatProductOperation(st->Op,MATPRODUCT_AB,NULL,MatMatMult_STOperator,NULL,MATDENSEHIP,MATDENSEHIP));
348:     }
349:     /* make sure the shell matrix generates a vector of the same type as the problem matrices */
350:     PetscCall(MatCreateVecs(st->A[0],&v,NULL));
351:     PetscCall(VecGetType(v,&vtype));
352:     PetscCall(MatShellSetVecType(st->Op,vtype));
353:     PetscCall(VecDestroy(&v));
354:     /* build the operator matrices */
355:     PetscCall(STComputeOperator(st));
356:   }
357:   if (Op) *Op = st->Op;
358:   PetscFunctionReturn(PETSC_SUCCESS);
359: }

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

365:    Collective

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

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

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

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

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

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

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

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

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

401:    Level: advanced

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

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

422:    Logically Collective

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

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

431:    Level: advanced

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

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

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

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

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

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

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

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

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

499:    Collective

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

504:    Level: advanced

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

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

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

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

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

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

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

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

662:    Collective

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

667:    Level: developer

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

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

685:    Not Collective

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

693:    Level: developer

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

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

711:    Not Collective

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

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

719:    Level: developer

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

727:   PetscFunctionBegin;
730:   PetscAssertPointer(is,2);

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

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

741:    Collective

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

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

756:    Level: developer

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

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

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

779:    Collective

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

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

788:    Level: developer

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

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