Actual source code: stfunc.c
 
   slepc-3.18.3 2023-03-24
   
  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:    The ST interface routines, callable by users
 12: */
 14: #include <slepc/private/stimpl.h>
 16: PetscClassId     ST_CLASSID = 0;
 17: PetscLogEvent    ST_SetUp = 0,ST_ComputeOperator = 0,ST_Apply = 0,ST_ApplyTranspose = 0,ST_MatSetUp = 0,ST_MatMult = 0,ST_MatMultTranspose = 0,ST_MatSolve = 0,ST_MatSolveTranspose = 0;
 18: static PetscBool STPackageInitialized = PETSC_FALSE;
 20: const char *STMatModes[] = {"COPY","INPLACE","SHELL","STMatMode","ST_MATMODE_",NULL};
 22: /*@C
 23:    STFinalizePackage - This function destroys everything in the Slepc interface
 24:    to the ST package. It is called from SlepcFinalize().
 26:    Level: developer
 28: .seealso: SlepcFinalize()
 29: @*/
 30: PetscErrorCode STFinalizePackage(void)
 31: {
 32:   PetscFunctionListDestroy(&STList);
 33:   STPackageInitialized = PETSC_FALSE;
 34:   STRegisterAllCalled  = PETSC_FALSE;
 35:   return 0;
 36: }
 38: /*@C
 39:    STInitializePackage - This function initializes everything in the ST package.
 40:    It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 41:    on the first call to STCreate() when using static libraries.
 43:    Level: developer
 45: .seealso: SlepcInitialize()
 46: @*/
 47: PetscErrorCode STInitializePackage(void)
 48: {
 49:   char           logList[256];
 50:   PetscBool      opt,pkg;
 51:   PetscClassId   classids[1];
 53:   if (STPackageInitialized) return 0;
 54:   STPackageInitialized = PETSC_TRUE;
 55:   /* Register Classes */
 56:   PetscClassIdRegister("Spectral Transform",&ST_CLASSID);
 57:   /* Register Constructors */
 58:   STRegisterAll();
 59:   /* Register Events */
 60:   PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp);
 61:   PetscLogEventRegister("STComputeOperatr",ST_CLASSID,&ST_ComputeOperator);
 62:   PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply);
 63:   PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose);
 64:   PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp);
 65:   PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult);
 66:   PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose);
 67:   PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve);
 68:   PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose);
 69:   /* Process Info */
 70:   classids[0] = ST_CLASSID;
 71:   PetscInfoProcessClass("st",1,&classids[0]);
 72:   /* Process summary exclusions */
 73:   PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);
 74:   if (opt) {
 75:     PetscStrInList("st",logList,',',&pkg);
 76:     if (pkg) PetscLogEventDeactivateClass(ST_CLASSID);
 77:   }
 78:   /* Register package finalizer */
 79:   PetscRegisterFinalize(STFinalizePackage);
 80:   return 0;
 81: }
 83: /*@
 84:    STReset - Resets the ST context to the initial state (prior to setup)
 85:    and destroys any allocated Vecs and Mats.
 87:    Collective on st
 89:    Input Parameter:
 90: .  st - the spectral transformation context
 92:    Level: advanced
 94: .seealso: STDestroy()
 95: @*/
 96: PetscErrorCode STReset(ST st)
 97: {
 99:   if (!st) return 0;
100:   STCheckNotSeized(st,1);
101:   PetscTryTypeMethod(st,reset);
102:   if (st->ksp) KSPReset(st->ksp);
103:   MatDestroyMatrices(PetscMax(2,st->nmat),&st->T);
104:   MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
105:   st->nmat = 0;
106:   PetscFree(st->Astate);
107:   MatDestroy(&st->Op);
108:   MatDestroy(&st->P);
109:   MatDestroy(&st->Pmat);
110:   MatDestroyMatrices(st->nsplit,&st->Psplit);
111:   st->nsplit = 0;
112:   VecDestroyVecs(st->nwork,&st->work);
113:   st->nwork = 0;
114:   VecDestroy(&st->wb);
115:   VecDestroy(&st->wht);
116:   VecDestroy(&st->D);
117:   st->state   = ST_STATE_INITIAL;
118:   st->opready = PETSC_FALSE;
119:   return 0;
120: }
122: /*@C
123:    STDestroy - Destroys ST context that was created with STCreate().
125:    Collective on st
127:    Input Parameter:
128: .  st - the spectral transformation context
130:    Level: beginner
132: .seealso: STCreate(), STSetUp()
133: @*/
134: PetscErrorCode STDestroy(ST *st)
135: {
136:   if (!*st) return 0;
138:   if (--((PetscObject)(*st))->refct > 0) { *st = NULL; return 0; }
139:   STReset(*st);
140:   PetscTryTypeMethod(*st,destroy);
141:   KSPDestroy(&(*st)->ksp);
142:   PetscHeaderDestroy(st);
143:   return 0;
144: }
146: /*@
147:    STCreate - Creates a spectral transformation context.
149:    Collective
151:    Input Parameter:
152: .  comm - MPI communicator
154:    Output Parameter:
155: .  newst - location to put the spectral transformation context
157:    Level: beginner
159: .seealso: STSetUp(), STApply(), STDestroy(), ST
160: @*/
161: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)
162: {
163:   ST             st;
166:   *newst = NULL;
167:   STInitializePackage();
168:   SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);
170:   st->A            = NULL;
171:   st->nmat         = 0;
172:   st->sigma        = 0.0;
173:   st->defsigma     = 0.0;
174:   st->matmode      = ST_MATMODE_COPY;
175:   st->str          = UNKNOWN_NONZERO_PATTERN;
176:   st->transform    = PETSC_FALSE;
177:   st->D            = NULL;
178:   st->Pmat         = NULL;
179:   st->Pmat_set     = PETSC_FALSE;
180:   st->Psplit       = NULL;
181:   st->nsplit       = 0;
182:   st->strp         = UNKNOWN_NONZERO_PATTERN;
184:   st->ksp          = NULL;
185:   st->usesksp      = PETSC_FALSE;
186:   st->nwork        = 0;
187:   st->work         = NULL;
188:   st->wb           = NULL;
189:   st->wht          = NULL;
190:   st->state        = ST_STATE_INITIAL;
191:   st->Astate       = NULL;
192:   st->T            = NULL;
193:   st->Op           = NULL;
194:   st->opseized     = PETSC_FALSE;
195:   st->opready      = PETSC_FALSE;
196:   st->P            = NULL;
197:   st->M            = NULL;
198:   st->sigma_set    = PETSC_FALSE;
199:   st->asymm        = PETSC_FALSE;
200:   st->aherm        = PETSC_FALSE;
201:   st->data         = NULL;
203:   *newst = st;
204:   return 0;
205: }
207: /*
208:    Checks whether the ST matrices are all symmetric or hermitian.
209: */
210: static inline PetscErrorCode STMatIsSymmetricKnown(ST st,PetscBool *symm,PetscBool *herm)
211: {
212:   PetscInt       i;
213:   PetscBool      sbaij=PETSC_FALSE,set,flg=PETSC_FALSE;
215:   /* check if problem matrices are all sbaij */
216:   for (i=0;i<st->nmat;i++) {
217:     PetscObjectTypeCompareAny((PetscObject)st->A[i],&sbaij,MATSEQSBAIJ,MATMPISBAIJ,"");
218:     if (!sbaij) break;
219:   }
220:   /* check if user has set the symmetric flag */
221:   *symm = PETSC_TRUE;
222:   for (i=0;i<st->nmat;i++) {
223:     MatIsSymmetricKnown(st->A[i],&set,&flg);
224:     if (!set || !flg) { *symm = PETSC_FALSE; break; }
225:   }
226:   if (sbaij) *symm = PETSC_TRUE;
227: #if defined(PETSC_USE_COMPLEX)
228:   /* check if user has set the hermitian flag */
229:   *herm = PETSC_TRUE;
230:   for (i=0;i<st->nmat;i++) {
231:     MatIsHermitianKnown(st->A[i],&set,&flg);
232:     if (!set || !flg) { *herm = PETSC_FALSE; break; }
233:   }
234: #else
235:   *herm = *symm;
236: #endif
237:   return 0;
238: }
240: /*@
241:    STSetMatrices - Sets the matrices associated with the eigenvalue problem.
243:    Collective on st
245:    Input Parameters:
246: +  st - the spectral transformation context
247: .  n  - number of matrices in array A
248: -  A  - the array of matrices associated with the eigensystem
250:    Notes:
251:    It must be called before STSetUp(). If it is called again after STSetUp() then
252:    the ST object is reset.
254:    Level: intermediate
256: .seealso: STGetMatrix(), STGetNumMatrices(), STSetUp(), STReset()
257: @*/
258: PetscErrorCode STSetMatrices(ST st,PetscInt n,Mat A[])
259: {
260:   PetscInt       i;
261:   PetscBool      same=PETSC_TRUE;
268:   STCheckNotSeized(st,1);
271:   if (st->state) {
272:     if (n!=st->nmat) same = PETSC_FALSE;
273:     for (i=0;same&&i<n;i++) {
274:       if (A[i]!=st->A[i]) same = PETSC_FALSE;
275:     }
276:     if (!same) STReset(st);
277:   } else same = PETSC_FALSE;
278:   if (!same) {
279:     MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
280:     PetscCalloc1(PetscMax(2,n),&st->A);
281:     PetscFree(st->Astate);
282:     PetscMalloc1(PetscMax(2,n),&st->Astate);
283:   }
284:   for (i=0;i<n;i++) {
286:     PetscObjectReference((PetscObject)A[i]);
287:     MatDestroy(&st->A[i]);
288:     st->A[i] = A[i];
289:     st->Astate[i] = ((PetscObject)A[i])->state;
290:   }
291:   if (n==1) {
292:     st->A[1] = NULL;
293:     st->Astate[1] = 0;
294:   }
295:   st->nmat = n;
296:   if (same) st->state = ST_STATE_UPDATED;
297:   else st->state = ST_STATE_INITIAL;
299:   st->opready = PETSC_FALSE;
300:   if (!same) STMatIsSymmetricKnown(st,&st->asymm,&st->aherm);
301:   return 0;
302: }
304: /*@
305:    STGetMatrix - Gets the matrices associated with the original eigensystem.
307:    Not collective, though parallel Mats are returned if the ST is parallel
309:    Input Parameters:
310: +  st - the spectral transformation context
311: -  k  - the index of the requested matrix (starting in 0)
313:    Output Parameters:
314: .  A - the requested matrix
316:    Level: intermediate
318: .seealso: STSetMatrices(), STGetNumMatrices()
319: @*/
320: PetscErrorCode STGetMatrix(ST st,PetscInt k,Mat *A)
321: {
325:   STCheckMatrices(st,1);
328:   *A = st->A[k];
329:   return 0;
330: }
332: /*@
333:    STGetMatrixTransformed - Gets the matrices associated with the transformed eigensystem.
335:    Not collective, though parallel Mats are returned if the ST is parallel
337:    Input Parameters:
338: +  st - the spectral transformation context
339: -  k  - the index of the requested matrix (starting in 0)
341:    Output Parameters:
342: .  T - the requested matrix
344:    Level: developer
346: .seealso: STGetMatrix(), STGetNumMatrices()
347: @*/
348: PetscErrorCode STGetMatrixTransformed(ST st,PetscInt k,Mat *T)
349: {
353:   STCheckMatrices(st,1);
356:   *T = st->T[k];
357:   return 0;
358: }
360: /*@
361:    STGetNumMatrices - Returns the number of matrices stored in the ST.
363:    Not collective
365:    Input Parameter:
366: .  st - the spectral transformation context
368:    Output Parameters:
369: .  n - the number of matrices passed in STSetMatrices()
371:    Level: intermediate
373: .seealso: STSetMatrices()
374: @*/
375: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)
376: {
379:   *n = st->nmat;
380:   return 0;
381: }
383: /*@
384:    STResetMatrixState - Resets the stored state of the matrices in the ST.
386:    Logically Collective on st
388:    Input Parameter:
389: .  st - the spectral transformation context
391:    Note:
392:    This is useful in solvers where the user matrices are modified during
393:    the computation, as in nonlinear inverse iteration. The effect is that
394:    STGetMatrix() will retrieve the modified matrices as if they were
395:    the matrices originally provided by the user.
397:    Level: developer
399: .seealso: STGetMatrix(), EPSPowerSetNonlinear()
400: @*/
401: PetscErrorCode STResetMatrixState(ST st)
402: {
403:   PetscInt i;
406:   for (i=0;i<st->nmat;i++) st->Astate[i] = ((PetscObject)st->A[i])->state;
407:   return 0;
408: }
410: /*@
411:    STSetPreconditionerMat - Sets the matrix to be used to build the preconditioner.
413:    Collective on st
415:    Input Parameters:
416: +  st  - the spectral transformation context
417: -  mat - the matrix that will be used in constructing the preconditioner
419:    Notes:
420:    This matrix will be passed to the internal KSP object (via the last argument
421:    of KSPSetOperators()) as the matrix to be used when constructing the preconditioner.
422:    If no matrix is set or mat is set to NULL, A-sigma*B will be used
423:    to build the preconditioner, being sigma the value set by STSetShift().
425:    More precisely, this is relevant for spectral transformations that represent
426:    a rational matrix function, and use a KSP object for the denominator, called
427:    K in the description of STGetOperator(). It includes also the STPRECOND case.
428:    If the user has a good approximation to matrix K that can be used to build a
429:    cheap preconditioner, it can be passed with this function. Note that it affects
430:    only the Pmat argument of KSPSetOperators(), not the Amat argument.
432:    If a preconditioner matrix is set, the default is to use an iterative KSP
433:    rather than a direct method.
435:    An alternative to pass an approximation of A-sigma*B with this function is
436:    to provide approximations of A and B via STSetSplitPreconditioner(). The
437:    difference is that when sigma changes the preconditioner is recomputed.
439:    Use NULL to remove a previously set matrix.
441:    Level: advanced
443: .seealso: STGetPreconditionerMat(), STSetShift(), STGetOperator(), STSetSplitPreconditioner()
444: @*/
445: PetscErrorCode STSetPreconditionerMat(ST st,Mat mat)
446: {
448:   if (mat) {
451:   }
452:   STCheckNotSeized(st,1);
454:   if (mat) PetscObjectReference((PetscObject)mat);
455:   MatDestroy(&st->Pmat);
456:   st->Pmat     = mat;
457:   st->Pmat_set = mat? PETSC_TRUE: PETSC_FALSE;
458:   st->state    = ST_STATE_INITIAL;
459:   st->opready  = PETSC_FALSE;
460:   return 0;
461: }
463: /*@
464:    STGetPreconditionerMat - Returns the matrix previously set by STSetPreconditionerMat().
466:    Not Collective
468:    Input Parameter:
469: .  st - the spectral transformation context
471:    Output Parameter:
472: .  mat - the matrix that will be used in constructing the preconditioner or
473:    NULL if no matrix was set by STSetPreconditionerMat().
475:    Level: advanced
477: .seealso: STSetPreconditionerMat()
478: @*/
479: PetscErrorCode STGetPreconditionerMat(ST st,Mat *mat)
480: {
483:   *mat = st->Pmat_set? st->Pmat: NULL;
484:   return 0;
485: }
487: /*@
488:    STSetSplitPreconditioner - Sets the matrices from which to build the preconditioner
489:    in split form.
491:    Collective on st
493:    Input Parameters:
494: +  st     - the spectral transformation context
495: .  n      - number of matrices
496: .  Psplit - array of matrices
497: -  strp   - structure flag for Psplit matrices
499:    Notes:
500:    The number of matrices passed here must be the same as in STSetMatrices().
502:    For linear eigenproblems, the preconditioner matrix is computed as
503:    Pmat(sigma) = A0-sigma*B0, where A0 and B0 are approximations of A and B
504:    (the eigenproblem matrices) provided via the Psplit array in this function.
505:    Compared to STSetPreconditionerMat(), this function allows setting a preconditioner
506:    in a way that is independent of the shift sigma. Whenever the value of sigma
507:    changes the preconditioner is recomputed.
509:    Similarly, for polynomial eigenproblems the matrix for the preconditioner
510:    is expressed as Pmat(sigma) = sum_i Psplit_i*phi_i(sigma), for i=1,...,n, where
511:    the phi_i's are the polynomial basis functions.
513:    The structure flag provides information about the relative nonzero pattern of the
514:    Psplit_i matrices, in the same way as in STSetMatStructure().
516:    Use n=0 to reset a previously set split preconditioner.
518:    Level: advanced
520: .seealso: STGetSplitPreconditionerTerm(), STGetSplitPreconditionerInfo(), STSetPreconditionerMat(), STSetMatrices(), STSetMatStructure()
521: @*/
522: PetscErrorCode STSetSplitPreconditioner(ST st,PetscInt n,Mat Psplit[],MatStructure strp)
523: {
524:   PetscInt       i,N=0,M,M0=0,mloc,nloc,mloc0=0;
533:   STCheckNotSeized(st,1);
535:   for (i=0;i<n;i++) {
538:     MatGetSize(Psplit[i],&M,&N);
539:     MatGetLocalSize(Psplit[i],&mloc,&nloc);
542:     if (!i) { M0 = M; mloc0 = mloc; }
545:     PetscObjectReference((PetscObject)Psplit[i]);
546:   }
548:   if (st->Psplit) MatDestroyMatrices(st->nsplit,&st->Psplit);
550:   /* allocate space and copy matrices */
551:   if (n) {
552:     PetscMalloc1(n,&st->Psplit);
553:     for (i=0;i<n;i++) st->Psplit[i] = Psplit[i];
554:   }
555:   st->nsplit = n;
556:   st->strp   = strp;
557:   st->state  = ST_STATE_INITIAL;
558:   return 0;
559: }
561: /*@
562:    STGetSplitPreconditionerTerm - Gets the matrices associated with
563:    the split preconditioner.
565:    Not collective, though parallel Mats are returned if the ST is parallel
567:    Input Parameters:
568: +  st - the spectral transformation context
569: -  k  - the index of the requested matrix (starting in 0)
571:    Output Parameter:
572: .  Psplit - the returned matrix
574:    Level: advanced
576: .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerInfo()
577: @*/
578: PetscErrorCode STGetSplitPreconditionerTerm(ST st,PetscInt k,Mat *Psplit)
579: {
585:   *Psplit = st->Psplit[k];
586:   return 0;
587: }
589: /*@
590:    STGetSplitPreconditionerInfo - Returns the number of matrices of the split
591:    preconditioner, as well as the structure flag.
593:    Not collective
595:    Input Parameter:
596: .  st - the spectral transformation context
598:    Output Parameters:
599: +  n    - the number of matrices passed in STSetSplitPreconditioner()
600: -  strp - the matrix structure flag passed in STSetSplitPreconditioner()
602:    Level: advanced
604: .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerTerm()
605: @*/
606: PetscErrorCode STGetSplitPreconditionerInfo(ST st,PetscInt *n,MatStructure *strp)
607: {
609:   if (n)    *n    = st->nsplit;
610:   if (strp) *strp = st->strp;
611:   return 0;
612: }
614: /*@
615:    STSetShift - Sets the shift associated with the spectral transformation.
617:    Logically Collective on st
619:    Input Parameters:
620: +  st - the spectral transformation context
621: -  shift - the value of the shift
623:    Notes:
624:    In some spectral transformations, changing the shift may have associated
625:    a lot of work, for example recomputing a factorization.
627:    This function is normally not directly called by users, since the shift is
628:    indirectly set by EPSSetTarget().
630:    Level: intermediate
632: .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift()
633: @*/
634: PetscErrorCode STSetShift(ST st,PetscScalar shift)
635: {
639:   if (st->sigma != shift) {
640:     STCheckNotSeized(st,1);
641:     if (st->state==ST_STATE_SETUP) PetscTryTypeMethod(st,setshift,shift);
642:     st->sigma = shift;
643:   }
644:   st->sigma_set = PETSC_TRUE;
645:   return 0;
646: }
648: /*@
649:    STGetShift - Gets the shift associated with the spectral transformation.
651:    Not Collective
653:    Input Parameter:
654: .  st - the spectral transformation context
656:    Output Parameter:
657: .  shift - the value of the shift
659:    Level: intermediate
661: .seealso: STSetShift()
662: @*/
663: PetscErrorCode STGetShift(ST st,PetscScalar* shift)
664: {
667:   *shift = st->sigma;
668:   return 0;
669: }
671: /*@
672:    STSetDefaultShift - Sets the value of the shift that should be employed if
673:    the user did not specify one.
675:    Logically Collective on st
677:    Input Parameters:
678: +  st - the spectral transformation context
679: -  defaultshift - the default value of the shift
681:    Level: developer
683: .seealso: STSetShift()
684: @*/
685: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)
686: {
689:   if (st->defsigma != defaultshift) {
690:     st->defsigma = defaultshift;
691:     st->state    = ST_STATE_INITIAL;
692:     st->opready  = PETSC_FALSE;
693:   }
694:   return 0;
695: }
697: /*@
698:    STScaleShift - Multiply the shift with a given factor.
700:    Logically Collective on st
702:    Input Parameters:
703: +  st     - the spectral transformation context
704: -  factor - the scaling factor
706:    Note:
707:    This function does not update the transformation matrices, as opposed to
708:    STSetShift().
710:    Level: developer
712: .seealso: STSetShift()
713: @*/
714: PetscErrorCode STScaleShift(ST st,PetscScalar factor)
715: {
718:   st->sigma *= factor;
719:   return 0;
720: }
722: /*@
723:    STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
725:    Collective on st
727:    Input Parameters:
728: +  st - the spectral transformation context
729: -  D  - the diagonal matrix (represented as a vector)
731:    Notes:
732:    If this matrix is set, STApply will effectively apply D*OP*D^{-1}. Use NULL
733:    to reset a previously passed D.
735:    Balancing is usually set via EPSSetBalance, but the advanced user may use
736:    this function to bypass the usual balancing methods.
738:    Level: developer
740: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
741: @*/
742: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)
743: {
745:   if (st->D == D) return 0;
746:   STCheckNotSeized(st,1);
747:   if (D) {
750:     PetscObjectReference((PetscObject)D);
751:   }
752:   VecDestroy(&st->D);
753:   st->D = D;
754:   st->state   = ST_STATE_INITIAL;
755:   st->opready = PETSC_FALSE;
756:   return 0;
757: }
759: /*@
760:    STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
762:    Not collective, but vector is shared by all processors that share the ST
764:    Input Parameter:
765: .  st - the spectral transformation context
767:    Output Parameter:
768: .  D  - the diagonal matrix (represented as a vector)
770:    Note:
771:    If the matrix was not set, a null pointer will be returned.
773:    Level: developer
775: .seealso: STSetBalanceMatrix()
776: @*/
777: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)
778: {
781:   *D = st->D;
782:   return 0;
783: }
785: /*@C
786:    STMatCreateVecs - Get vector(s) compatible with the ST matrices.
788:    Collective on st
790:    Input Parameter:
791: .  st - the spectral transformation context
793:    Output Parameters:
794: +  right - (optional) vector that the matrix can be multiplied against
795: -  left  - (optional) vector that the matrix vector product can be stored in
797:    Level: developer
799: .seealso: STMatCreateVecsEmpty()
800: @*/
801: PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left)
802: {
803:   STCheckMatrices(st,1);
804:   MatCreateVecs(st->A[0],right,left);
805:   return 0;
806: }
808: /*@C
809:    STMatCreateVecsEmpty - Get vector(s) compatible with the ST matrices, i.e. with the same
810:    parallel layout, but without internal array.
812:    Collective on st
814:    Input Parameter:
815: .  st - the spectral transformation context
817:    Output Parameters:
818: +  right - (optional) vector that the matrix can be multiplied against
819: -  left  - (optional) vector that the matrix vector product can be stored in
821:    Level: developer
823: .seealso: STMatCreateVecs(), MatCreateVecsEmpty()
824: @*/
825: PetscErrorCode STMatCreateVecsEmpty(ST st,Vec *right,Vec *left)
826: {
827:   STCheckMatrices(st,1);
828:   MatCreateVecsEmpty(st->A[0],right,left);
829:   return 0;
830: }
832: /*@
833:    STMatGetSize - Returns the number of rows and columns of the ST matrices.
835:    Not Collective
837:    Input Parameter:
838: .  st - the spectral transformation context
840:    Output Parameters:
841: +  m - the number of global rows
842: -  n - the number of global columns
844:    Level: developer
846: .seealso: STMatGetLocalSize()
847: @*/
848: PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)
849: {
850:   STCheckMatrices(st,1);
851:   MatGetSize(st->A[0],m,n);
852:   return 0;
853: }
855: /*@
856:    STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.
858:    Not Collective
860:    Input Parameter:
861: .  st - the spectral transformation context
863:    Output Parameters:
864: +  m - the number of local rows
865: -  n - the number of local columns
867:    Level: developer
869: .seealso: STMatGetSize()
870: @*/
871: PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)
872: {
873:   STCheckMatrices(st,1);
874:   MatGetLocalSize(st->A[0],m,n);
875:   return 0;
876: }
878: /*@C
879:    STSetOptionsPrefix - Sets the prefix used for searching for all
880:    ST options in the database.
882:    Logically Collective on st
884:    Input Parameters:
885: +  st     - the spectral transformation context
886: -  prefix - the prefix string to prepend to all ST option requests
888:    Notes:
889:    A hyphen (-) must NOT be given at the beginning of the prefix name.
890:    The first character of all runtime options is AUTOMATICALLY the
891:    hyphen.
893:    Level: advanced
895: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
896: @*/
897: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)
898: {
900:   if (!st->ksp) STGetKSP(st,&st->ksp);
901:   KSPSetOptionsPrefix(st->ksp,prefix);
902:   KSPAppendOptionsPrefix(st->ksp,"st_");
903:   PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
904:   return 0;
905: }
907: /*@C
908:    STAppendOptionsPrefix - Appends to the prefix used for searching for all
909:    ST options in the database.
911:    Logically Collective on st
913:    Input Parameters:
914: +  st     - the spectral transformation context
915: -  prefix - the prefix string to prepend to all ST option requests
917:    Notes:
918:    A hyphen (-) must NOT be given at the beginning of the prefix name.
919:    The first character of all runtime options is AUTOMATICALLY the
920:    hyphen.
922:    Level: advanced
924: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
925: @*/
926: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)
927: {
929:   PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
930:   if (!st->ksp) STGetKSP(st,&st->ksp);
931:   KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
932:   KSPAppendOptionsPrefix(st->ksp,"st_");
933:   return 0;
934: }
936: /*@C
937:    STGetOptionsPrefix - Gets the prefix used for searching for all
938:    ST options in the database.
940:    Not Collective
942:    Input Parameters:
943: .  st - the spectral transformation context
945:    Output Parameters:
946: .  prefix - pointer to the prefix string used, is returned
948:    Note:
949:    On the Fortran side, the user should pass in a string 'prefix' of
950:    sufficient length to hold the prefix.
952:    Level: advanced
954: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
955: @*/
956: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])
957: {
960:   PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
961:   return 0;
962: }
964: /*@C
965:    STView - Prints the ST data structure.
967:    Collective on st
969:    Input Parameters:
970: +  st - the ST context
971: -  viewer - optional visualization context
973:    Note:
974:    The available visualization contexts include
975: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
976: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
977:          output where only the first processor opens
978:          the file.  All other processors send their
979:          data to the first processor to print.
981:    The user can open an alternative visualization contexts with
982:    PetscViewerASCIIOpen() (output to a specified file).
984:    Level: beginner
986: .seealso: EPSView()
987: @*/
988: PetscErrorCode STView(ST st,PetscViewer viewer)
989: {
990:   STType         cstr;
991:   char           str[50];
992:   PetscBool      isascii,isstring;
995:   if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)st),&viewer);
999:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1000:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1001:   if (isascii) {
1002:     PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer);
1003:     PetscViewerASCIIPushTab(viewer);
1004:     PetscTryTypeMethod(st,view,viewer);
1005:     PetscViewerASCIIPopTab(viewer);
1006:     SlepcSNPrintfScalar(str,sizeof(str),st->sigma,PETSC_FALSE);
1007:     PetscViewerASCIIPrintf(viewer,"  shift: %s\n",str);
1008:     PetscViewerASCIIPrintf(viewer,"  number of matrices: %" PetscInt_FMT "\n",st->nmat);
1009:     switch (st->matmode) {
1010:     case ST_MATMODE_COPY:
1011:       break;
1012:     case ST_MATMODE_INPLACE:
1013:       PetscViewerASCIIPrintf(viewer,"  shifting the matrix and unshifting at exit\n");
1014:       break;
1015:     case ST_MATMODE_SHELL:
1016:       PetscViewerASCIIPrintf(viewer,"  using a shell matrix\n");
1017:       break;
1018:     }
1019:     if (st->nmat>1 && st->matmode != ST_MATMODE_SHELL) PetscViewerASCIIPrintf(viewer,"  nonzero pattern of the matrices: %s\n",MatStructures[st->str]);
1020:     if (st->Psplit) PetscViewerASCIIPrintf(viewer,"  using split preconditioner matrices with %s\n",MatStructures[st->strp]);
1021:     if (st->transform && st->nmat>2) PetscViewerASCIIPrintf(viewer,"  computing transformed matrices\n");
1022:   } else if (isstring) {
1023:     STGetType(st,&cstr);
1024:     PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1025:     PetscTryTypeMethod(st,view,viewer);
1026:   }
1027:   if (st->usesksp) {
1028:     if (!st->ksp) STGetKSP(st,&st->ksp);
1029:     PetscViewerASCIIPushTab(viewer);
1030:     KSPView(st->ksp,viewer);
1031:     PetscViewerASCIIPopTab(viewer);
1032:   }
1033:   return 0;
1034: }
1036: /*@C
1037:    STViewFromOptions - View from options
1039:    Collective on ST
1041:    Input Parameters:
1042: +  st   - the spectral transformation context
1043: .  obj  - optional object
1044: -  name - command line option
1046:    Level: intermediate
1048: .seealso: STView(), STCreate()
1049: @*/
1050: PetscErrorCode STViewFromOptions(ST st,PetscObject obj,const char name[])
1051: {
1053:   PetscObjectViewFromOptions((PetscObject)st,obj,name);
1054:   return 0;
1055: }
1057: /*@C
1058:    STRegister - Adds a method to the spectral transformation package.
1060:    Not collective
1062:    Input Parameters:
1063: +  name - name of a new user-defined transformation
1064: -  function - routine to create method context
1066:    Notes:
1067:    STRegister() may be called multiple times to add several user-defined
1068:    spectral transformations.
1070:    Sample usage:
1071: .vb
1072:     STRegister("my_transform",MyTransformCreate);
1073: .ve
1075:    Then, your spectral transform can be chosen with the procedural interface via
1076: $     STSetType(st,"my_transform")
1077:    or at runtime via the option
1078: $     -st_type my_transform
1080:    Level: advanced
1082: .seealso: STRegisterAll()
1083: @*/
1084: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))
1085: {
1086:   STInitializePackage();
1087:   PetscFunctionListAdd(&STList,name,function);
1088:   return 0;
1089: }