Actual source code: stsolve.c

slepc-3.13.0 2020-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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>            /*I "slepcst.h" I*/

 16: PetscErrorCode STApply_Generic(ST st,Vec x,Vec y)
 17: {

 21:   if (st->M && st->P) {
 22:     MatMult(st->M,x,st->work[0]);
 23:     STMatSolve(st,st->work[0],y);
 24:   } else if (st->M) {
 25:     MatMult(st->M,x,y);
 26:   } else {
 27:     STMatSolve(st,x,y);
 28:   }
 29:   return(0);
 30: }

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

 37:    Collective on st

 39:    Input Parameters:
 40: +  st - the spectral transformation context
 41: -  x  - input vector

 43:    Output Parameter:
 44: .  y - output vector

 46:    Level: developer

 48: .seealso: STApplyTranspose(), STApplyHermitianTranspose()
 49: @*/
 50: PetscErrorCode STApply(ST st,Vec x,Vec y)
 51: {
 53:   Mat            Op;

 60:   STCheckMatrices(st,1);
 61:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
 62:   VecSetErrorIfLocked(y,3);
 63:   if (!st->ops->apply) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST does not have apply");
 64:   STGetOperator_Private(st,&Op);
 65:   MatMult(Op,x,y);
 66:   return(0);
 67: }

 69: PetscErrorCode STApplyTranspose_Generic(ST st,Vec x,Vec y)
 70: {

 74:   if (st->M && st->P) {
 75:     STMatSolveTranspose(st,x,st->work[0]);
 76:     MatMultTranspose(st->M,st->work[0],y);
 77:   } else if (st->M) {
 78:     MatMultTranspose(st->M,x,y);
 79:   } else {
 80:     STMatSolveTranspose(st,x,y);
 81:   }
 82:   return(0);
 83: }

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

 90:    Collective on st

 92:    Input Parameters:
 93: +  st - the spectral transformation context
 94: -  x  - input vector

 96:    Output Parameter:
 97: .  y - output vector

 99:    Level: developer

101: .seealso: STApply(), STApplyHermitianTranspose()
102: @*/
103: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
104: {
106:   Mat            Op;

113:   STCheckMatrices(st,1);
114:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
115:   VecSetErrorIfLocked(y,3);
116:   if (!st->ops->applytrans) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST does not have applytrans");
117:   STGetOperator_Private(st,&Op);
118:   MatMultTranspose(Op,x,y);
119:   return(0);
120: }

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

127:    Collective on st

129:    Input Parameters:
130: +  st - the spectral transformation context
131: -  x  - input vector

133:    Output Parameter:
134: .  y - output vector

136:    Note:
137:    Currently implemented via STApplyTranspose() with appropriate conjugation.

139:    Level: developer

141: .seealso: STApply(), STApplyTranspose()
142: @*/
143: PetscErrorCode STApplyHermitianTranspose(ST st,Vec x,Vec y)
144: {
146:   Mat            Op;

153:   STCheckMatrices(st,1);
154:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
155:   VecSetErrorIfLocked(y,3);
156:   if (!st->ops->applytrans) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST does not have applytrans");
157:   STGetOperator_Private(st,&Op);
158:   MatMultHermitianTranspose(Op,x,y);
159:   return(0);
160: }

162: /*@
163:    STGetBilinearForm - Returns the matrix used in the bilinear form with a
164:    generalized problem with semi-definite B.

166:    Not collective, though a parallel Mat may be returned

168:    Input Parameters:
169: .  st - the spectral transformation context

171:    Output Parameter:
172: .  B - output matrix

174:    Notes:
175:    The output matrix B must be destroyed after use. It will be NULL in
176:    case of standard eigenproblems.

178:    Level: developer
179: @*/
180: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
181: {

188:   STCheckMatrices(st,1);
189:   (*st->ops->getbilinearform)(st,B);
190:   return(0);
191: }

193: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
194: {

198:   if (st->nmat==1) *B = NULL;
199:   else {
200:     *B = st->A[1];
201:     PetscObjectReference((PetscObject)*B);
202:   }
203:   return(0);
204: }

206: static PetscErrorCode MatMult_STOperator(Mat Op,Vec x,Vec y)
207: {
209:   ST             st;

212:   MatShellGetContext(Op,(void**)&st);
213:   STSetUp(st);
214:   PetscLogEventBegin(ST_Apply,st,x,y,0);
215:   if (st->D) { /* with balancing */
216:     VecPointwiseDivide(st->wb,x,st->D);
217:     (*st->ops->apply)(st,st->wb,y);
218:     VecPointwiseMult(y,y,st->D);
219:   } else {
220:     (*st->ops->apply)(st,x,y);
221:   }
222:   PetscLogEventEnd(ST_Apply,st,x,y,0);
223:   return(0);
224: }

226: static PetscErrorCode MatMultTranspose_STOperator(Mat Op,Vec x,Vec y)
227: {
229:   ST             st;

232:   MatShellGetContext(Op,(void**)&st);
233:   STSetUp(st);
234:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
235:   if (st->D) { /* with balancing */
236:     VecPointwiseMult(st->wb,x,st->D);
237:     (*st->ops->applytrans)(st,st->wb,y);
238:     VecPointwiseDivide(y,y,st->D);
239:   } else {
240:     (*st->ops->applytrans)(st,x,y);
241:   }
242:   PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
243:   return(0);
244: }

246: #if defined(PETSC_USE_COMPLEX)
247: static PetscErrorCode MatMultHermitianTranspose_STOperator(Mat Op,Vec x,Vec y)
248: {
250:   ST             st;

253:   MatShellGetContext(Op,(void**)&st);
254:   STSetUp(st);
255:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
256:   if (!st->wht) {
257:     MatCreateVecs(st->A[0],&st->wht,NULL);
258:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->wht);
259:   }
260:   VecCopy(x,st->wht);
261:   VecConjugate(st->wht);
262:   if (st->D) { /* with balancing */
263:     VecPointwiseMult(st->wb,st->wht,st->D);
264:     (*st->ops->applytrans)(st,st->wb,y);
265:     VecPointwiseDivide(y,y,st->D);
266:   } else {
267:     (*st->ops->applytrans)(st,st->wht,y);
268:   }
269:   VecConjugate(y);
270:   PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
271:   return(0);
272: }
273: #endif

275: PetscErrorCode STGetOperator_Private(ST st,Mat *Op)
276: {
278:   PetscInt       m,n,M,N;

280:   if (!st->Op) {
281:     MatGetLocalSize(st->A[0],&m,&n);
282:     MatGetSize(st->A[0],&M,&N);
283:     MatCreateShell(PetscObjectComm((PetscObject)st),m,n,M,N,st,&st->Op);
284:     MatShellSetOperation(st->Op,MATOP_MULT,(void(*)(void))MatMult_STOperator);
285:     MatShellSetOperation(st->Op,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_STOperator);
286: #if defined(PETSC_USE_COMPLEX)
287:     MatShellSetOperation(st->Op,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultHermitianTranspose_STOperator);
288: #else
289:     MatShellSetOperation(st->Op,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultTranspose_STOperator);
290:     STComputeOperator(st);
291: #endif
292:   }
293:   if (Op) *Op = st->Op;
294:   return(0);
295: }

297: /*@
298:    STGetOperator - Returns a shell matrix that represents the operator of the
299:    spectral transformation.

301:    Collective on st

303:    Input Parameter:
304: .  st - the spectral transformation context

306:    Output Parameter:
307: .  Op - operator matrix

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

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

316: .vb
317:       Op = D*inv(K)*M*inv(D)
318: .ve

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

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

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

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

339:    Level: advanced

341: .seealso: STApply(), STApplyTranspose(), STSetBalanceMatrix(), STShellSetApply(),
342:           STGetKSP(), STSetShift(), STRestoreOperator(), STSetMatrices()
343: @*/
344: PetscErrorCode STGetOperator(ST st,Mat *Op)
345: {

351:   STCheckMatrices(st,1);
352:   STCheckNotSeized(st,1);
353:   if (st->nmat>2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"The operator is not defined in polynomial eigenproblems");
354:   STGetOperator_Private(st,Op);
355:   if (Op) st->opseized = PETSC_TRUE;
356:   return(0);
357: }

359: /*@
360:    STRestoreOperator - Restore the previously seized operator matrix.

362:    Collective on st

364:    Input Parameters:
365: +  st - the spectral transformation context
366: -  Op - operator matrix

368:    Notes:
369:    The arguments must match the corresponding call to STGetOperator().

371:    Level: advanced

373: .seealso: STGetOperator()
374: @*/
375: PetscErrorCode STRestoreOperator(ST st,Mat *Op)
376: {
381:   if (!st->opseized) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Must be called after STGetOperator()");
382:   *Op = NULL;
383:   st->opseized = PETSC_FALSE;
384:   return(0);
385: }

387: /*
388:    STComputeOperator - Computes the matrices that constitute the operator

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

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

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

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

404:    STSHELL does not compute anything here, but sets the flag as if it was ready.
405: */
406: PetscErrorCode STComputeOperator(ST st)
407: {
409:   PC             pc;

414:   if (!st->opready && st->ops->computeoperator) {
415:     STCheckMatrices(st,1);
416:     if (!st->T) {
417:       PetscCalloc1(PetscMax(2,st->nmat),&st->T);
418:       PetscLogObjectMemory((PetscObject)st,PetscMax(2,st->nmat)*sizeof(Mat));
419:     }
420:     PetscLogEventBegin(ST_ComputeOperator,st,0,0,0);
421:     (*st->ops->computeoperator)(st);
422:     PetscLogEventEnd(ST_ComputeOperator,st,0,0,0);
423:     if (st->usesksp) {
424:       if (!st->ksp) { STGetKSP(st,&st->ksp); }
425:       if (st->P) {
426:         STSetDefaultKSP(st);
427:         STCheckFactorPackage(st);
428:         KSPSetOperators(st->ksp,st->P,st->P);
429:       } else {
430:         /* STPRECOND defaults to PCNONE if st->P is empty */
431:         KSPGetPC(st->ksp,&pc);
432:         PCSetType(pc,PCNONE);
433:       }
434:     }
435:   }
436:   st->opready = PETSC_TRUE;
437:   return(0);
438: }

440: /*@
441:    STSetUp - Prepares for the use of a spectral transformation.

443:    Collective on st

445:    Input Parameter:
446: .  st - the spectral transformation context

448:    Level: advanced

450: .seealso: STCreate(), STApply(), STDestroy()
451: @*/
452: PetscErrorCode STSetUp(ST st)
453: {
454:   PetscInt       i,n,k;

460:   STCheckMatrices(st,1);
461:   switch (st->state) {
462:     case ST_STATE_INITIAL:
463:       PetscInfo(st,"Setting up new ST\n");
464:       if (!((PetscObject)st)->type_name) {
465:         STSetType(st,STSHIFT);
466:       }
467:       break;
468:     case ST_STATE_SETUP:
469:       return(0);
470:     case ST_STATE_UPDATED:
471:       PetscInfo(st,"Setting up updated ST\n");
472:       break;
473:   }
474:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
475:   if (st->state!=ST_STATE_UPDATED) {
476:     if (!(st->nmat<3 && st->opready)) {
477:       if (st->T) {
478:         for (i=0;i<PetscMax(2,st->nmat);i++) {
479:           MatDestroy(&st->T[i]);
480:         }
481:       }
482:       MatDestroy(&st->P);
483:     }
484:   }
485:   if (st->D) {
486:     MatGetLocalSize(st->A[0],NULL,&n);
487:     VecGetLocalSize(st->D,&k);
488:     if (n != k) SETERRQ2(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_SIZ,"Balance matrix has wrong dimension %D (should be %D)",k,n);
489:     if (!st->wb) {
490:       VecDuplicate(st->D,&st->wb);
491:       PetscLogObjectParent((PetscObject)st,(PetscObject)st->wb);
492:     }
493:   }
494:   if (st->nmat<3 && st->transform) {
495:     STComputeOperator(st);
496:   } else {
497:     if (!st->T) {
498:       PetscCalloc1(PetscMax(2,st->nmat),&st->T);
499:       PetscLogObjectMemory((PetscObject)st,PetscMax(2,st->nmat)*sizeof(Mat));
500:     }
501:   }
502:   if (st->ops->setup) { (*st->ops->setup)(st); }
503:   st->state = ST_STATE_SETUP;
504:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
505:   return(0);
506: }

508: /*
509:    Computes coefficients for the transformed polynomial,
510:    and stores the result in argument S.

512:    alpha - value of the parameter of the transformed polynomial
513:    beta - value of the previous shift (only used in inplace mode)
514:    k - index of first matrix included in the computation
515:    coeffs - coefficients of the expansion
516:    initial - true if this is the first time (only relevant for shell mode)
517: */
518: PetscErrorCode STMatMAXPY_Private(ST st,PetscScalar alpha,PetscScalar beta,PetscInt k,PetscScalar *coeffs,PetscBool initial,Mat *S)
519: {
521:   PetscInt       *matIdx=NULL,nmat,i,ini=-1;
522:   PetscScalar    t=1.0,ta,gamma;
523:   PetscBool      nz=PETSC_FALSE;

526:   nmat = st->nmat-k;
527:   switch (st->matmode) {
528:   case ST_MATMODE_INPLACE:
529:     if (st->nmat>2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST_MATMODE_INPLACE not supported for polynomial eigenproblems");
530:     if (initial) {
531:       PetscObjectReference((PetscObject)st->A[0]);
532:       *S = st->A[0];
533:       gamma = alpha;
534:     } else gamma = alpha-beta;
535:     if (gamma != 0.0) {
536:       if (st->nmat>1) {
537:         MatAXPY(*S,gamma,st->A[1],st->str);
538:       } else {
539:         MatShift(*S,gamma);
540:       }
541:     }
542:     break;
543:   case ST_MATMODE_SHELL:
544:     if (initial) {
545:       if (st->nmat>2) {
546:         PetscMalloc1(nmat,&matIdx);
547:         for (i=0;i<nmat;i++) matIdx[i] = k+i;
548:       }
549:       STMatShellCreate(st,alpha,nmat,matIdx,coeffs,S);
550:       PetscLogObjectParent((PetscObject)st,(PetscObject)*S);
551:       if (st->nmat>2) { PetscFree(matIdx); }
552:     } else {
553:       STMatShellShift(*S,alpha);
554:     }
555:     break;
556:   case ST_MATMODE_COPY:
557:     if (coeffs) {
558:       for (i=0;i<nmat && ini==-1;i++) {
559:         if (coeffs[i]!=0.0) ini = i;
560:         else t *= alpha;
561:       }
562:       if (coeffs[ini] != 1.0) nz = PETSC_TRUE;
563:       for (i=ini+1;i<nmat&&!nz;i++) if (coeffs[i]!=0.0) nz = PETSC_TRUE;
564:     } else { nz = PETSC_TRUE; ini = 0; }
565:     if ((alpha == 0.0 || !nz) && t==1.0) {
566:       PetscObjectReference((PetscObject)st->A[k+ini]);
567:       MatDestroy(S);
568:       *S = st->A[k+ini];
569:     } else {
570:       if (*S && *S!=st->A[k+ini]) {
571:         MatSetOption(*S,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
572:         MatCopy(st->A[k+ini],*S,DIFFERENT_NONZERO_PATTERN);
573:       } else {
574:         MatDestroy(S);
575:         MatDuplicate(st->A[k+ini],MAT_COPY_VALUES,S);
576:         MatSetOption(*S,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
577:         PetscLogObjectParent((PetscObject)st,(PetscObject)*S);
578:       }
579:       if (coeffs && coeffs[ini]!=1.0) {
580:         MatScale(*S,coeffs[ini]);
581:       }
582:       for (i=ini+k+1;i<PetscMax(2,st->nmat);i++) {
583:         t *= alpha;
584:         ta = t;
585:         if (coeffs) ta *= coeffs[i-k];
586:         if (ta!=0.0) {
587:           if (st->nmat>1) {
588:             MatAXPY(*S,ta,st->A[i],st->str);
589:           } else {
590:             MatShift(*S,ta);
591:           }
592:         }
593:       }
594:     }
595:   }
596:   MatSetOption(*S,MAT_SYMMETRIC,st->asymm);
597:   MatSetOption(*S,MAT_HERMITIAN,(PetscImaginaryPart(st->sigma)==0.0)?st->aherm:PETSC_FALSE);
598:   return(0);
599: }

601: /*
602:    Computes the values of the coefficients required by STMatMAXPY_Private
603:    for the case of monomial basis.
604: */
605: PetscErrorCode STCoeffs_Monomial(ST st, PetscScalar *coeffs)
606: {
607:   PetscInt  k,i,ini,inip;

610:   /* Compute binomial coefficients */
611:   ini = (st->nmat*(st->nmat-1))/2;
612:   for (i=0;i<st->nmat;i++) coeffs[ini+i]=1.0;
613:   for (k=st->nmat-1;k>=1;k--) {
614:     inip = ini+1;
615:     ini = (k*(k-1))/2;
616:     coeffs[ini] = 1.0;
617:     for (i=1;i<k;i++) coeffs[ini+i] = coeffs[ini+i-1]+coeffs[inip+i-1];
618:   }
619:   return(0);
620: }

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

626:    Collective on st

628:    Input Parameters:
629: .  st  - the spectral transformation context

631:    Level: developer

633: .seealso: EPSSolve()
634: @*/
635: PetscErrorCode STPostSolve(ST st)
636: {

642:   if (st->ops->postsolve) {
643:     (*st->ops->postsolve)(st);
644:   }
645:   return(0);
646: }

648: /*@
649:    STBackTransform - Back-transformation phase, intended for
650:    spectral transformations which require to transform the computed
651:    eigenvalues back to the original eigenvalue problem.

653:    Not Collective

655:    Input Parameters:
656: +  st   - the spectral transformation context
657:    eigr - real part of a computed eigenvalue
658: -  eigi - imaginary part of a computed eigenvalue

660:    Level: developer

662: .seealso: STIsInjective()
663: @*/
664: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
665: {

671:   if (st->ops->backtransform) {
672:     (*st->ops->backtransform)(st,n,eigr,eigi);
673:   }
674:   return(0);
675: }

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

682:    Not collective

684:    Input Parameter:
685: .  st   - the spectral transformation context

687:    Output Parameter:
688: .  is - the answer

690:    Level: developer

692: .seealso: STBackTransform()
693: @*/
694: PetscErrorCode STIsInjective(ST st,PetscBool* is)
695: {
697:   PetscBool      shell;


704:   PetscObjectTypeCompare((PetscObject)st,STSHELL,&shell);
705:   if (shell) {
706:     STIsInjective_Shell(st,is);
707:   } else *is = st->ops->backtransform? PETSC_TRUE: PETSC_FALSE;
708:   return(0);
709: }

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

714:    Collective on st

716:    Input Parameters:
717: +  st     - the spectral transformation context
718: .  sigma  - the shift
719: -  coeffs - the coefficients (may be NULL)

721:    Note:
722:    This function is not intended to be called by end users, but by SLEPc
723:    solvers that use ST. It builds matrix st->P as follows, then calls KSPSetUp().
724: .vb
725:     If (coeffs):  st->P = Sum_{i=0:nmat-1} coeffs[i]*sigma^i*A_i.
726:     else          st->P = Sum_{i=0:nmat-1} sigma^i*A_i
727: .ve

729:    Level: developer

731: .seealso: STMatSolve()
732: @*/
733: PetscErrorCode STMatSetUp(ST st,PetscScalar sigma,PetscScalar *coeffs)
734: {

740:   STCheckMatrices(st,1);

742:   PetscLogEventBegin(ST_MatSetUp,st,0,0,0);
743:   STMatMAXPY_Private(st,sigma,0.0,0,coeffs,PETSC_TRUE,&st->P);
744:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
745:   STCheckFactorPackage(st);
746:   KSPSetOperators(st->ksp,st->P,st->P);
747:   KSPSetUp(st->ksp);
748:   PetscLogEventEnd(ST_MatSetUp,st,0,0,0);
749:   return(0);
750: }

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

755:    Collective on st

757:    Input Parameters:
758: +  st - the spectral transformation context
759: -  nw - number of work vectors to allocate

761:    Developers Note:
762:    This is SLEPC_EXTERN because it may be required by shell STs.

764:    Level: developer
765: @*/
766: PetscErrorCode STSetWorkVecs(ST st,PetscInt nw)
767: {
769:   PetscInt       i;

774:   if (nw <= 0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"nw must be > 0: nw = %D",nw);
775:   if (st->nwork < nw) {
776:     VecDestroyVecs(st->nwork,&st->work);
777:     st->nwork = nw;
778:     PetscMalloc1(nw,&st->work);
779:     for (i=0;i<nw;i++) { STMatCreateVecs(st,&st->work[i],NULL); }
780:     PetscLogObjectParents(st,nw,st->work);
781:   }
782:   return(0);
783: }