LCOV - code coverage report
Current view: top level - sys/classes/st/interface - stfunc.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 375 399 94.0 %
Date: 2024-04-18 00:38:59 Functions: 29 32 90.6 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       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             : */
      13             : 
      14             : #include <slepc/private/stimpl.h>            /*I "slepcst.h" I*/
      15             : 
      16             : PetscClassId     ST_CLASSID = 0;
      17             : PetscLogEvent    ST_SetUp = 0,ST_ComputeOperator = 0,ST_Apply = 0,ST_ApplyTranspose = 0,ST_ApplyHermitianTranspose = 0,ST_MatSetUp = 0,ST_MatMult = 0,ST_MatMultTranspose = 0,ST_MatSolve = 0,ST_MatSolveTranspose = 0;
      18             : static PetscBool STPackageInitialized = PETSC_FALSE;
      19             : 
      20             : const char *STMatModes[] = {"COPY","INPLACE","SHELL","STMatMode","ST_MATMODE_",NULL};
      21             : 
      22             : /*@C
      23             :    STFinalizePackage - This function destroys everything in the Slepc interface
      24             :    to the ST package. It is called from SlepcFinalize().
      25             : 
      26             :    Level: developer
      27             : 
      28             : .seealso: SlepcFinalize()
      29             : @*/
      30         723 : PetscErrorCode STFinalizePackage(void)
      31             : {
      32         723 :   PetscFunctionBegin;
      33         723 :   PetscCall(PetscFunctionListDestroy(&STList));
      34         723 :   STPackageInitialized = PETSC_FALSE;
      35         723 :   STRegisterAllCalled  = PETSC_FALSE;
      36         723 :   PetscFunctionReturn(PETSC_SUCCESS);
      37             : }
      38             : 
      39             : /*@C
      40             :    STInitializePackage - This function initializes everything in the ST package.
      41             :    It is called from PetscDLLibraryRegister() when using dynamic libraries, and
      42             :    on the first call to STCreate() when using static libraries.
      43             : 
      44             :    Level: developer
      45             : 
      46             : .seealso: SlepcInitialize()
      47             : @*/
      48        5147 : PetscErrorCode STInitializePackage(void)
      49             : {
      50        5147 :   char           logList[256];
      51        5147 :   PetscBool      opt,pkg;
      52        5147 :   PetscClassId   classids[1];
      53             : 
      54        5147 :   PetscFunctionBegin;
      55        5147 :   if (STPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
      56         723 :   STPackageInitialized = PETSC_TRUE;
      57             :   /* Register Classes */
      58         723 :   PetscCall(PetscClassIdRegister("Spectral Transform",&ST_CLASSID));
      59             :   /* Register Constructors */
      60         723 :   PetscCall(STRegisterAll());
      61             :   /* Register Events */
      62         723 :   PetscCall(PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp));
      63         723 :   PetscCall(PetscLogEventRegister("STComputeOperatr",ST_CLASSID,&ST_ComputeOperator));
      64         723 :   PetscCall(PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply));
      65         723 :   PetscCall(PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose));
      66         723 :   PetscCall(PetscLogEventRegister("STApplyHermTrans",ST_CLASSID,&ST_ApplyHermitianTranspose));
      67         723 :   PetscCall(PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp));
      68         723 :   PetscCall(PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult));
      69         723 :   PetscCall(PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose));
      70         723 :   PetscCall(PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve));
      71         723 :   PetscCall(PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose));
      72             :   /* Process Info */
      73         723 :   classids[0] = ST_CLASSID;
      74         723 :   PetscCall(PetscInfoProcessClass("st",1,&classids[0]));
      75             :   /* Process summary exclusions */
      76         723 :   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
      77         723 :   if (opt) {
      78           8 :     PetscCall(PetscStrInList("st",logList,',',&pkg));
      79           8 :     if (pkg) PetscCall(PetscLogEventDeactivateClass(ST_CLASSID));
      80             :   }
      81             :   /* Register package finalizer */
      82         723 :   PetscCall(PetscRegisterFinalize(STFinalizePackage));
      83         723 :   PetscFunctionReturn(PETSC_SUCCESS);
      84             : }
      85             : 
      86             : /*@
      87             :    STReset - Resets the ST context to the initial state (prior to setup)
      88             :    and destroys any allocated Vecs and Mats.
      89             : 
      90             :    Collective
      91             : 
      92             :    Input Parameter:
      93             : .  st - the spectral transformation context
      94             : 
      95             :    Level: advanced
      96             : 
      97             : .seealso: STDestroy()
      98             : @*/
      99        1810 : PetscErrorCode STReset(ST st)
     100             : {
     101        1810 :   PetscFunctionBegin;
     102        1810 :   if (st) PetscValidHeaderSpecific(st,ST_CLASSID,1);
     103           0 :   if (!st) PetscFunctionReturn(PETSC_SUCCESS);
     104        1810 :   STCheckNotSeized(st,1);
     105        1810 :   PetscTryTypeMethod(st,reset);
     106        1810 :   if (st->ksp) PetscCall(KSPReset(st->ksp));
     107        1810 :   PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->T));
     108        1810 :   PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->A));
     109        1810 :   st->nmat = 0;
     110        1810 :   PetscCall(PetscFree(st->Astate));
     111        1810 :   PetscCall(MatDestroy(&st->Op));
     112        1810 :   PetscCall(MatDestroy(&st->P));
     113        1810 :   PetscCall(MatDestroy(&st->Pmat));
     114        1810 :   PetscCall(MatDestroyMatrices(st->nsplit,&st->Psplit));
     115        1810 :   st->nsplit = 0;
     116        1810 :   PetscCall(VecDestroyVecs(st->nwork,&st->work));
     117        1810 :   st->nwork = 0;
     118        1810 :   PetscCall(VecDestroy(&st->wb));
     119        1810 :   PetscCall(VecDestroy(&st->wht));
     120        1810 :   PetscCall(VecDestroy(&st->D));
     121        1810 :   st->state   = ST_STATE_INITIAL;
     122        1810 :   st->opready = PETSC_FALSE;
     123        1810 :   PetscFunctionReturn(PETSC_SUCCESS);
     124             : }
     125             : 
     126             : /*@C
     127             :    STDestroy - Destroys ST context that was created with STCreate().
     128             : 
     129             :    Collective
     130             : 
     131             :    Input Parameter:
     132             : .  st - the spectral transformation context
     133             : 
     134             :    Level: beginner
     135             : 
     136             : .seealso: STCreate(), STSetUp()
     137             : @*/
     138         812 : PetscErrorCode STDestroy(ST *st)
     139             : {
     140         812 :   PetscFunctionBegin;
     141         812 :   if (!*st) PetscFunctionReturn(PETSC_SUCCESS);
     142         810 :   PetscValidHeaderSpecific(*st,ST_CLASSID,1);
     143         810 :   if (--((PetscObject)*st)->refct > 0) { *st = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
     144         808 :   PetscCall(STReset(*st));
     145         808 :   PetscTryTypeMethod(*st,destroy);
     146         808 :   PetscCall(KSPDestroy(&(*st)->ksp));
     147         808 :   PetscCall(PetscHeaderDestroy(st));
     148         808 :   PetscFunctionReturn(PETSC_SUCCESS);
     149             : }
     150             : 
     151             : /*@
     152             :    STCreate - Creates a spectral transformation context.
     153             : 
     154             :    Collective
     155             : 
     156             :    Input Parameter:
     157             : .  comm - MPI communicator
     158             : 
     159             :    Output Parameter:
     160             : .  newst - location to put the spectral transformation context
     161             : 
     162             :    Level: beginner
     163             : 
     164             : .seealso: STSetUp(), STApply(), STDestroy(), ST
     165             : @*/
     166         808 : PetscErrorCode STCreate(MPI_Comm comm,ST *newst)
     167             : {
     168         808 :   ST             st;
     169             : 
     170         808 :   PetscFunctionBegin;
     171         808 :   PetscAssertPointer(newst,2);
     172         808 :   *newst = NULL;
     173         808 :   PetscCall(STInitializePackage());
     174         808 :   PetscCall(SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView));
     175             : 
     176         808 :   st->A            = NULL;
     177         808 :   st->nmat         = 0;
     178         808 :   st->sigma        = 0.0;
     179         808 :   st->defsigma     = 0.0;
     180         808 :   st->matmode      = ST_MATMODE_COPY;
     181         808 :   st->str          = UNKNOWN_NONZERO_PATTERN;
     182         808 :   st->transform    = PETSC_FALSE;
     183         808 :   st->D            = NULL;
     184         808 :   st->Pmat         = NULL;
     185         808 :   st->Pmat_set     = PETSC_FALSE;
     186         808 :   st->Psplit       = NULL;
     187         808 :   st->nsplit       = 0;
     188         808 :   st->strp         = UNKNOWN_NONZERO_PATTERN;
     189             : 
     190         808 :   st->ksp          = NULL;
     191         808 :   st->usesksp      = PETSC_FALSE;
     192         808 :   st->nwork        = 0;
     193         808 :   st->work         = NULL;
     194         808 :   st->wb           = NULL;
     195         808 :   st->wht          = NULL;
     196         808 :   st->state        = ST_STATE_INITIAL;
     197         808 :   st->Astate       = NULL;
     198         808 :   st->T            = NULL;
     199         808 :   st->Op           = NULL;
     200         808 :   st->opseized     = PETSC_FALSE;
     201         808 :   st->opready      = PETSC_FALSE;
     202         808 :   st->P            = NULL;
     203         808 :   st->M            = NULL;
     204         808 :   st->sigma_set    = PETSC_FALSE;
     205         808 :   st->asymm        = PETSC_FALSE;
     206         808 :   st->aherm        = PETSC_FALSE;
     207         808 :   st->data         = NULL;
     208             : 
     209         808 :   *newst = st;
     210         808 :   PetscFunctionReturn(PETSC_SUCCESS);
     211             : }
     212             : 
     213             : /*
     214             :    Checks whether the ST matrices are all symmetric or hermitian.
     215             : */
     216         889 : static inline PetscErrorCode STMatIsSymmetricKnown(ST st,PetscBool *symm,PetscBool *herm)
     217             : {
     218         889 :   PetscInt       i;
     219         889 :   PetscBool      sbaij=PETSC_FALSE,set,flg=PETSC_FALSE;
     220             : 
     221         889 :   PetscFunctionBegin;
     222             :   /* check if problem matrices are all sbaij */
     223         893 :   for (i=0;i<st->nmat;i++) {
     224         891 :     PetscCall(PetscObjectTypeCompareAny((PetscObject)st->A[i],&sbaij,MATSEQSBAIJ,MATMPISBAIJ,""));
     225         891 :     if (!sbaij) break;
     226             :   }
     227             :   /* check if user has set the symmetric flag */
     228         889 :   *symm = PETSC_TRUE;
     229         949 :   for (i=0;i<st->nmat;i++) {
     230         919 :     PetscCall(MatIsSymmetricKnown(st->A[i],&set,&flg));
     231         919 :     if (!set || !flg) { *symm = PETSC_FALSE; break; }
     232             :   }
     233         889 :   if (sbaij) *symm = PETSC_TRUE;
     234             : #if defined(PETSC_USE_COMPLEX)
     235             :   /* check if user has set the hermitian flag */
     236         889 :   *herm = PETSC_TRUE;
     237         954 :   for (i=0;i<st->nmat;i++) {
     238         917 :     PetscCall(MatIsHermitianKnown(st->A[i],&set,&flg));
     239         917 :     if (!set || !flg) { *herm = PETSC_FALSE; break; }
     240             :   }
     241             : #else
     242             :   *herm = *symm;
     243             : #endif
     244         889 :   PetscFunctionReturn(PETSC_SUCCESS);
     245             : }
     246             : 
     247             : /*@
     248             :    STSetMatrices - Sets the matrices associated with the eigenvalue problem.
     249             : 
     250             :    Collective
     251             : 
     252             :    Input Parameters:
     253             : +  st - the spectral transformation context
     254             : .  n  - number of matrices in array A
     255             : -  A  - the array of matrices associated with the eigensystem
     256             : 
     257             :    Notes:
     258             :    It must be called before STSetUp(). If it is called again after STSetUp() then
     259             :    the ST object is reset.
     260             : 
     261             :    Level: intermediate
     262             : 
     263             : .seealso: STGetMatrix(), STGetNumMatrices(), STSetUp(), STReset()
     264             : @*/
     265         930 : PetscErrorCode STSetMatrices(ST st,PetscInt n,Mat A[])
     266             : {
     267         930 :   PetscInt       i;
     268         930 :   PetscBool      same=PETSC_TRUE;
     269             : 
     270         930 :   PetscFunctionBegin;
     271         930 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     272        3720 :   PetscValidLogicalCollectiveInt(st,n,2);
     273         930 :   PetscCheck(n>0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %" PetscInt_FMT,n);
     274         930 :   PetscAssertPointer(A,3);
     275         930 :   PetscCheckSameComm(st,1,*A,3);
     276         930 :   STCheckNotSeized(st,1);
     277         930 :   PetscCheck(!st->nsplit || st->nsplit==n,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"The number of matrices must be the same as in STSetSplitPreconditioner()");
     278             : 
     279         930 :   if (st->state) {
     280          76 :     if (n!=st->nmat) same = PETSC_FALSE;
     281         165 :     for (i=0;same&&i<n;i++) {
     282          89 :       if (A[i]!=st->A[i]) same = PETSC_FALSE;
     283             :     }
     284          76 :     if (!same) PetscCall(STReset(st));
     285             :   } else same = PETSC_FALSE;
     286             :   if (!same) {
     287         889 :     PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->A));
     288         889 :     PetscCall(PetscCalloc1(PetscMax(2,n),&st->A));
     289         889 :     PetscCall(PetscFree(st->Astate));
     290         889 :     PetscCall(PetscMalloc1(PetscMax(2,n),&st->Astate));
     291             :   }
     292        2630 :   for (i=0;i<n;i++) {
     293        1700 :     PetscValidHeaderSpecific(A[i],MAT_CLASSID,3);
     294        1700 :     PetscCall(PetscObjectReference((PetscObject)A[i]));
     295        1700 :     PetscCall(MatDestroy(&st->A[i]));
     296        1700 :     st->A[i] = A[i];
     297        1700 :     st->Astate[i] = ((PetscObject)A[i])->state;
     298             :   }
     299         930 :   if (n==1) {
     300         522 :     st->A[1] = NULL;
     301         522 :     st->Astate[1] = 0;
     302             :   }
     303         930 :   st->nmat = n;
     304         930 :   if (same) st->state = ST_STATE_UPDATED;
     305         889 :   else st->state = ST_STATE_INITIAL;
     306         930 :   PetscCheck(!same || !st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Support for changing the matrices while using a split preconditioner is not implemented yet");
     307         930 :   st->opready = PETSC_FALSE;
     308         930 :   if (!same) PetscCall(STMatIsSymmetricKnown(st,&st->asymm,&st->aherm));
     309         930 :   PetscFunctionReturn(PETSC_SUCCESS);
     310             : }
     311             : 
     312             : /*@
     313             :    STGetMatrix - Gets the matrices associated with the original eigensystem.
     314             : 
     315             :    Not Collective
     316             : 
     317             :    Input Parameters:
     318             : +  st - the spectral transformation context
     319             : -  k  - the index of the requested matrix (starting in 0)
     320             : 
     321             :    Output Parameters:
     322             : .  A - the requested matrix
     323             : 
     324             :    Level: intermediate
     325             : 
     326             : .seealso: STSetMatrices(), STGetNumMatrices()
     327             : @*/
     328        9669 : PetscErrorCode STGetMatrix(ST st,PetscInt k,Mat *A)
     329             : {
     330        9669 :   PetscFunctionBegin;
     331        9669 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     332       38676 :   PetscValidLogicalCollectiveInt(st,k,2);
     333        9669 :   PetscAssertPointer(A,3);
     334        9669 :   STCheckMatrices(st,1);
     335        9669 :   PetscCheck(k>=0 && k<st->nmat,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat-1);
     336        9669 :   PetscCheck(((PetscObject)st->A[k])->state==st->Astate[k],PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
     337        9669 :   *A = st->A[k];
     338        9669 :   PetscFunctionReturn(PETSC_SUCCESS);
     339             : }
     340             : 
     341             : /*@
     342             :    STGetMatrixTransformed - Gets the matrices associated with the transformed eigensystem.
     343             : 
     344             :    Not Collective
     345             : 
     346             :    Input Parameters:
     347             : +  st - the spectral transformation context
     348             : -  k  - the index of the requested matrix (starting in 0)
     349             : 
     350             :    Output Parameters:
     351             : .  T - the requested matrix
     352             : 
     353             :    Level: developer
     354             : 
     355             : .seealso: STGetMatrix(), STGetNumMatrices()
     356             : @*/
     357         211 : PetscErrorCode STGetMatrixTransformed(ST st,PetscInt k,Mat *T)
     358             : {
     359         211 :   PetscFunctionBegin;
     360         211 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     361         844 :   PetscValidLogicalCollectiveInt(st,k,2);
     362         211 :   PetscAssertPointer(T,3);
     363         211 :   STCheckMatrices(st,1);
     364         211 :   PetscCheck(k>=0 && k<st->nmat,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat-1);
     365         211 :   PetscCheck(st->T,PetscObjectComm((PetscObject)st),PETSC_ERR_POINTER,"There are no transformed matrices");
     366         211 :   *T = st->T[k];
     367         211 :   PetscFunctionReturn(PETSC_SUCCESS);
     368             : }
     369             : 
     370             : /*@
     371             :    STGetNumMatrices - Returns the number of matrices stored in the ST.
     372             : 
     373             :    Not Collective
     374             : 
     375             :    Input Parameter:
     376             : .  st - the spectral transformation context
     377             : 
     378             :    Output Parameters:
     379             : .  n - the number of matrices passed in STSetMatrices()
     380             : 
     381             :    Level: intermediate
     382             : 
     383             : .seealso: STSetMatrices()
     384             : @*/
     385        7103 : PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)
     386             : {
     387        7103 :   PetscFunctionBegin;
     388        7103 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     389        7103 :   PetscAssertPointer(n,2);
     390        7103 :   *n = st->nmat;
     391        7103 :   PetscFunctionReturn(PETSC_SUCCESS);
     392             : }
     393             : 
     394             : /*@
     395             :    STResetMatrixState - Resets the stored state of the matrices in the ST.
     396             : 
     397             :    Logically Collective
     398             : 
     399             :    Input Parameter:
     400             : .  st - the spectral transformation context
     401             : 
     402             :    Note:
     403             :    This is useful in solvers where the user matrices are modified during
     404             :    the computation, as in nonlinear inverse iteration. The effect is that
     405             :    STGetMatrix() will retrieve the modified matrices as if they were
     406             :    the matrices originally provided by the user.
     407             : 
     408             :    Level: developer
     409             : 
     410             : .seealso: STGetMatrix(), EPSPowerSetNonlinear()
     411             : @*/
     412        1788 : PetscErrorCode STResetMatrixState(ST st)
     413             : {
     414        1788 :   PetscInt i;
     415             : 
     416        1788 :   PetscFunctionBegin;
     417        1788 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     418        5364 :   for (i=0;i<st->nmat;i++) st->Astate[i] = ((PetscObject)st->A[i])->state;
     419        1788 :   PetscFunctionReturn(PETSC_SUCCESS);
     420             : }
     421             : 
     422             : /*@
     423             :    STSetPreconditionerMat - Sets the matrix to be used to build the preconditioner.
     424             : 
     425             :    Collective
     426             : 
     427             :    Input Parameters:
     428             : +  st  - the spectral transformation context
     429             : -  mat - the matrix that will be used in constructing the preconditioner
     430             : 
     431             :    Notes:
     432             :    This matrix will be passed to the internal KSP object (via the last argument
     433             :    of KSPSetOperators()) as the matrix to be used when constructing the preconditioner.
     434             :    If no matrix is set or mat is set to NULL, A-sigma*B will be used
     435             :    to build the preconditioner, being sigma the value set by STSetShift().
     436             : 
     437             :    More precisely, this is relevant for spectral transformations that represent
     438             :    a rational matrix function, and use a KSP object for the denominator, called
     439             :    K in the description of STGetOperator(). It includes also the STPRECOND case.
     440             :    If the user has a good approximation to matrix K that can be used to build a
     441             :    cheap preconditioner, it can be passed with this function. Note that it affects
     442             :    only the Pmat argument of KSPSetOperators(), not the Amat argument.
     443             : 
     444             :    If a preconditioner matrix is set, the default is to use an iterative KSP
     445             :    rather than a direct method.
     446             : 
     447             :    An alternative to pass an approximation of A-sigma*B with this function is
     448             :    to provide approximations of A and B via STSetSplitPreconditioner(). The
     449             :    difference is that when sigma changes the preconditioner is recomputed.
     450             : 
     451             :    Use NULL to remove a previously set matrix.
     452             : 
     453             :    Level: advanced
     454             : 
     455             : .seealso: STGetPreconditionerMat(), STSetShift(), STGetOperator(), STSetSplitPreconditioner()
     456             : @*/
     457          21 : PetscErrorCode STSetPreconditionerMat(ST st,Mat mat)
     458             : {
     459          21 :   PetscFunctionBegin;
     460          21 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     461          21 :   if (mat) {
     462          20 :     PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
     463          20 :     PetscCheckSameComm(st,1,mat,2);
     464             :   }
     465          21 :   STCheckNotSeized(st,1);
     466          21 :   PetscCheck(!mat || !st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot call both STSetPreconditionerMat and STSetSplitPreconditioner");
     467          21 :   if (mat) PetscCall(PetscObjectReference((PetscObject)mat));
     468          21 :   PetscCall(MatDestroy(&st->Pmat));
     469          21 :   st->Pmat     = mat;
     470          21 :   st->Pmat_set = mat? PETSC_TRUE: PETSC_FALSE;
     471          21 :   st->state    = ST_STATE_INITIAL;
     472          21 :   st->opready  = PETSC_FALSE;
     473          21 :   PetscFunctionReturn(PETSC_SUCCESS);
     474             : }
     475             : 
     476             : /*@
     477             :    STGetPreconditionerMat - Returns the matrix previously set by STSetPreconditionerMat().
     478             : 
     479             :    Not Collective
     480             : 
     481             :    Input Parameter:
     482             : .  st - the spectral transformation context
     483             : 
     484             :    Output Parameter:
     485             : .  mat - the matrix that will be used in constructing the preconditioner or
     486             :    NULL if no matrix was set by STSetPreconditionerMat().
     487             : 
     488             :    Level: advanced
     489             : 
     490             : .seealso: STSetPreconditionerMat()
     491             : @*/
     492          43 : PetscErrorCode STGetPreconditionerMat(ST st,Mat *mat)
     493             : {
     494          43 :   PetscFunctionBegin;
     495          43 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     496          43 :   PetscAssertPointer(mat,2);
     497          43 :   *mat = st->Pmat_set? st->Pmat: NULL;
     498          43 :   PetscFunctionReturn(PETSC_SUCCESS);
     499             : }
     500             : 
     501             : /*@
     502             :    STSetSplitPreconditioner - Sets the matrices from which to build the preconditioner
     503             :    in split form.
     504             : 
     505             :    Collective
     506             : 
     507             :    Input Parameters:
     508             : +  st     - the spectral transformation context
     509             : .  n      - number of matrices
     510             : .  Psplit - array of matrices
     511             : -  strp   - structure flag for Psplit matrices
     512             : 
     513             :    Notes:
     514             :    The number of matrices passed here must be the same as in STSetMatrices().
     515             : 
     516             :    For linear eigenproblems, the preconditioner matrix is computed as
     517             :    Pmat(sigma) = A0-sigma*B0, where A0 and B0 are approximations of A and B
     518             :    (the eigenproblem matrices) provided via the Psplit array in this function.
     519             :    Compared to STSetPreconditionerMat(), this function allows setting a preconditioner
     520             :    in a way that is independent of the shift sigma. Whenever the value of sigma
     521             :    changes the preconditioner is recomputed.
     522             : 
     523             :    Similarly, for polynomial eigenproblems the matrix for the preconditioner
     524             :    is expressed as Pmat(sigma) = sum_i Psplit_i*phi_i(sigma), for i=1,...,n, where
     525             :    the phi_i's are the polynomial basis functions.
     526             : 
     527             :    The structure flag provides information about the relative nonzero pattern of the
     528             :    Psplit_i matrices, in the same way as in STSetMatStructure().
     529             : 
     530             :    Use n=0 to reset a previously set split preconditioner.
     531             : 
     532             :    Level: advanced
     533             : 
     534             : .seealso: STGetSplitPreconditionerTerm(), STGetSplitPreconditionerInfo(), STSetPreconditionerMat(), STSetMatrices(), STSetMatStructure()
     535             : @*/
     536          23 : PetscErrorCode STSetSplitPreconditioner(ST st,PetscInt n,Mat Psplit[],MatStructure strp)
     537             : {
     538          23 :   PetscInt       i,N=0,M,M0=0,mloc,nloc,mloc0=0;
     539             : 
     540          23 :   PetscFunctionBegin;
     541          23 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     542          92 :   PetscValidLogicalCollectiveInt(st,n,2);
     543          23 :   PetscCheck(n>=0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Negative value of n = %" PetscInt_FMT,n);
     544          23 :   PetscCheck(!n || !st->Pmat_set,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot call both STSetPreconditionerMat and STSetSplitPreconditioner");
     545          23 :   PetscCheck(!n || !st->nmat || st->nmat==n,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"The number of matrices must be the same as in STSetMatrices()");
     546          23 :   if (n) PetscAssertPointer(Psplit,3);
     547          92 :   PetscValidLogicalCollectiveEnum(st,strp,4);
     548          23 :   STCheckNotSeized(st,1);
     549             : 
     550          76 :   for (i=0;i<n;i++) {
     551          53 :     PetscValidHeaderSpecific(Psplit[i],MAT_CLASSID,3);
     552          53 :     PetscCheckSameComm(st,1,Psplit[i],3);
     553          53 :     PetscCall(MatGetSize(Psplit[i],&M,&N));
     554          53 :     PetscCall(MatGetLocalSize(Psplit[i],&mloc,&nloc));
     555          53 :     PetscCheck(M==N,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Psplit[%" PetscInt_FMT "] is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",i,M,N);
     556          53 :     PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Psplit[%" PetscInt_FMT "] does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,nloc);
     557          53 :     if (!i) { M0 = M; mloc0 = mloc; }
     558          53 :     PetscCheck(M==M0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_INCOMP,"Dimensions of Psplit[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,M,M0);
     559          53 :     PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_INCOMP,"Local dimensions of Psplit[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,mloc0);
     560          53 :     PetscCall(PetscObjectReference((PetscObject)Psplit[i]));
     561             :   }
     562             : 
     563          23 :   if (st->Psplit) PetscCall(MatDestroyMatrices(st->nsplit,&st->Psplit));
     564             : 
     565             :   /* allocate space and copy matrices */
     566          23 :   if (n) {
     567          23 :     PetscCall(PetscMalloc1(n,&st->Psplit));
     568          76 :     for (i=0;i<n;i++) st->Psplit[i] = Psplit[i];
     569             :   }
     570          23 :   st->nsplit = n;
     571          23 :   st->strp   = strp;
     572          23 :   st->state  = ST_STATE_INITIAL;
     573          23 :   PetscFunctionReturn(PETSC_SUCCESS);
     574             : }
     575             : 
     576             : /*@
     577             :    STGetSplitPreconditionerTerm - Gets the matrices associated with
     578             :    the split preconditioner.
     579             : 
     580             :    Not Collective
     581             : 
     582             :    Input Parameters:
     583             : +  st - the spectral transformation context
     584             : -  k  - the index of the requested matrix (starting in 0)
     585             : 
     586             :    Output Parameter:
     587             : .  Psplit - the returned matrix
     588             : 
     589             :    Level: advanced
     590             : 
     591             : .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerInfo()
     592             : @*/
     593          19 : PetscErrorCode STGetSplitPreconditionerTerm(ST st,PetscInt k,Mat *Psplit)
     594             : {
     595          19 :   PetscFunctionBegin;
     596          19 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     597          76 :   PetscValidLogicalCollectiveInt(st,k,2);
     598          19 :   PetscAssertPointer(Psplit,3);
     599          19 :   PetscCheck(k>=0 && k<st->nsplit,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nsplit-1);
     600          19 :   PetscCheck(st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_ORDER,"You have not called STSetSplitPreconditioner()");
     601          19 :   *Psplit = st->Psplit[k];
     602          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     603             : }
     604             : 
     605             : /*@
     606             :    STGetSplitPreconditionerInfo - Returns the number of matrices of the split
     607             :    preconditioner, as well as the structure flag.
     608             : 
     609             :    Not Collective
     610             : 
     611             :    Input Parameter:
     612             : .  st - the spectral transformation context
     613             : 
     614             :    Output Parameters:
     615             : +  n    - the number of matrices passed in STSetSplitPreconditioner()
     616             : -  strp - the matrix structure flag passed in STSetSplitPreconditioner()
     617             : 
     618             :    Level: advanced
     619             : 
     620             : .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerTerm()
     621             : @*/
     622         205 : PetscErrorCode STGetSplitPreconditionerInfo(ST st,PetscInt *n,MatStructure *strp)
     623             : {
     624         205 :   PetscFunctionBegin;
     625         205 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     626         205 :   if (n)    *n    = st->nsplit;
     627         205 :   if (strp) *strp = st->strp;
     628         205 :   PetscFunctionReturn(PETSC_SUCCESS);
     629             : }
     630             : 
     631             : /*@
     632             :    STSetShift - Sets the shift associated with the spectral transformation.
     633             : 
     634             :    Collective
     635             : 
     636             :    Input Parameters:
     637             : +  st - the spectral transformation context
     638             : -  shift - the value of the shift
     639             : 
     640             :    Notes:
     641             :    In some spectral transformations, changing the shift may have associated
     642             :    a lot of work, for example recomputing a factorization.
     643             : 
     644             :    This function is normally not directly called by users, since the shift is
     645             :    indirectly set by EPSSetTarget().
     646             : 
     647             :    Level: intermediate
     648             : 
     649             : .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift()
     650             : @*/
     651         716 : PetscErrorCode STSetShift(ST st,PetscScalar shift)
     652             : {
     653         716 :   PetscFunctionBegin;
     654         716 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     655         716 :   PetscValidType(st,1);
     656        2864 :   PetscValidLogicalCollectiveScalar(st,shift,2);
     657         716 :   if (st->sigma != shift) {
     658         698 :     STCheckNotSeized(st,1);
     659         698 :     if (st->state==ST_STATE_SETUP) PetscTryTypeMethod(st,setshift,shift);
     660         698 :     st->sigma = shift;
     661             :   }
     662         716 :   st->sigma_set = PETSC_TRUE;
     663         716 :   PetscFunctionReturn(PETSC_SUCCESS);
     664             : }
     665             : 
     666             : /*@
     667             :    STGetShift - Gets the shift associated with the spectral transformation.
     668             : 
     669             :    Not Collective
     670             : 
     671             :    Input Parameter:
     672             : .  st - the spectral transformation context
     673             : 
     674             :    Output Parameter:
     675             : .  shift - the value of the shift
     676             : 
     677             :    Level: intermediate
     678             : 
     679             : .seealso: STSetShift()
     680             : @*/
     681        3166 : PetscErrorCode STGetShift(ST st,PetscScalar* shift)
     682             : {
     683        3166 :   PetscFunctionBegin;
     684        3166 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     685        3166 :   PetscAssertPointer(shift,2);
     686        3166 :   *shift = st->sigma;
     687        3166 :   PetscFunctionReturn(PETSC_SUCCESS);
     688             : }
     689             : 
     690             : /*@
     691             :    STSetDefaultShift - Sets the value of the shift that should be employed if
     692             :    the user did not specify one.
     693             : 
     694             :    Logically Collective
     695             : 
     696             :    Input Parameters:
     697             : +  st - the spectral transformation context
     698             : -  defaultshift - the default value of the shift
     699             : 
     700             :    Level: developer
     701             : 
     702             : .seealso: STSetShift()
     703             : @*/
     704         306 : PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)
     705             : {
     706         306 :   PetscFunctionBegin;
     707         306 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     708        1224 :   PetscValidLogicalCollectiveScalar(st,defaultshift,2);
     709         306 :   if (st->defsigma != defaultshift) {
     710         202 :     st->defsigma = defaultshift;
     711         202 :     st->state    = ST_STATE_INITIAL;
     712         202 :     st->opready  = PETSC_FALSE;
     713             :   }
     714         306 :   PetscFunctionReturn(PETSC_SUCCESS);
     715             : }
     716             : 
     717             : /*@
     718             :    STScaleShift - Multiply the shift with a given factor.
     719             : 
     720             :    Logically Collective
     721             : 
     722             :    Input Parameters:
     723             : +  st     - the spectral transformation context
     724             : -  factor - the scaling factor
     725             : 
     726             :    Note:
     727             :    This function does not update the transformation matrices, as opposed to
     728             :    STSetShift().
     729             : 
     730             :    Level: developer
     731             : 
     732             : .seealso: STSetShift()
     733             : @*/
     734         216 : PetscErrorCode STScaleShift(ST st,PetscScalar factor)
     735             : {
     736         216 :   PetscFunctionBegin;
     737         216 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     738         864 :   PetscValidLogicalCollectiveScalar(st,factor,2);
     739         216 :   st->sigma *= factor;
     740         216 :   PetscFunctionReturn(PETSC_SUCCESS);
     741             : }
     742             : 
     743             : /*@
     744             :    STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
     745             : 
     746             :    Collective
     747             : 
     748             :    Input Parameters:
     749             : +  st - the spectral transformation context
     750             : -  D  - the diagonal matrix (represented as a vector)
     751             : 
     752             :    Notes:
     753             :    If this matrix is set, STApply will effectively apply D*OP*D^{-1}. Use NULL
     754             :    to reset a previously passed D.
     755             : 
     756             :    Balancing is usually set via EPSSetBalance, but the advanced user may use
     757             :    this function to bypass the usual balancing methods.
     758             : 
     759             :    Level: developer
     760             : 
     761             : .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
     762             : @*/
     763         861 : PetscErrorCode STSetBalanceMatrix(ST st,Vec D)
     764             : {
     765         861 :   PetscFunctionBegin;
     766         861 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     767         861 :   if (st->D == D) PetscFunctionReturn(PETSC_SUCCESS);
     768          13 :   STCheckNotSeized(st,1);
     769          13 :   if (D) {
     770          13 :     PetscValidHeaderSpecific(D,VEC_CLASSID,2);
     771          13 :     PetscCheckSameComm(st,1,D,2);
     772          13 :     PetscCall(PetscObjectReference((PetscObject)D));
     773             :   }
     774          13 :   PetscCall(VecDestroy(&st->D));
     775          13 :   st->D = D;
     776          13 :   st->state   = ST_STATE_INITIAL;
     777          13 :   st->opready = PETSC_FALSE;
     778          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     779             : }
     780             : 
     781             : /*@
     782             :    STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
     783             : 
     784             :    Not Collective
     785             : 
     786             :    Input Parameter:
     787             : .  st - the spectral transformation context
     788             : 
     789             :    Output Parameter:
     790             : .  D  - the diagonal matrix (represented as a vector)
     791             : 
     792             :    Note:
     793             :    If the matrix was not set, a null pointer will be returned.
     794             : 
     795             :    Level: developer
     796             : 
     797             : .seealso: STSetBalanceMatrix()
     798             : @*/
     799           0 : PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)
     800             : {
     801           0 :   PetscFunctionBegin;
     802           0 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     803           0 :   PetscAssertPointer(D,2);
     804           0 :   *D = st->D;
     805           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     806             : }
     807             : 
     808             : /*@C
     809             :    STMatCreateVecs - Get vector(s) compatible with the ST matrices.
     810             : 
     811             :    Collective
     812             : 
     813             :    Input Parameter:
     814             : .  st - the spectral transformation context
     815             : 
     816             :    Output Parameters:
     817             : +  right - (optional) vector that the matrix can be multiplied against
     818             : -  left  - (optional) vector that the matrix vector product can be stored in
     819             : 
     820             :    Level: developer
     821             : 
     822             : .seealso: STMatCreateVecsEmpty()
     823             : @*/
     824         380 : PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left)
     825             : {
     826         380 :   PetscFunctionBegin;
     827         380 :   STCheckMatrices(st,1);
     828         380 :   PetscCall(MatCreateVecs(st->A[0],right,left));
     829         380 :   PetscFunctionReturn(PETSC_SUCCESS);
     830             : }
     831             : 
     832             : /*@C
     833             :    STMatCreateVecsEmpty - Get vector(s) compatible with the ST matrices, i.e. with the same
     834             :    parallel layout, but without internal array.
     835             : 
     836             :    Collective
     837             : 
     838             :    Input Parameter:
     839             : .  st - the spectral transformation context
     840             : 
     841             :    Output Parameters:
     842             : +  right - (optional) vector that the matrix can be multiplied against
     843             : -  left  - (optional) vector that the matrix vector product can be stored in
     844             : 
     845             :    Level: developer
     846             : 
     847             : .seealso: STMatCreateVecs(), MatCreateVecsEmpty()
     848             : @*/
     849         826 : PetscErrorCode STMatCreateVecsEmpty(ST st,Vec *right,Vec *left)
     850             : {
     851         826 :   PetscFunctionBegin;
     852         826 :   STCheckMatrices(st,1);
     853         826 :   PetscCall(MatCreateVecsEmpty(st->A[0],right,left));
     854         826 :   PetscFunctionReturn(PETSC_SUCCESS);
     855             : }
     856             : 
     857             : /*@
     858             :    STMatGetSize - Returns the number of rows and columns of the ST matrices.
     859             : 
     860             :    Not Collective
     861             : 
     862             :    Input Parameter:
     863             : .  st - the spectral transformation context
     864             : 
     865             :    Output Parameters:
     866             : +  m - the number of global rows
     867             : -  n - the number of global columns
     868             : 
     869             :    Level: developer
     870             : 
     871             : .seealso: STMatGetLocalSize()
     872             : @*/
     873         848 : PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)
     874             : {
     875         848 :   PetscFunctionBegin;
     876         848 :   STCheckMatrices(st,1);
     877         848 :   PetscCall(MatGetSize(st->A[0],m,n));
     878         848 :   PetscFunctionReturn(PETSC_SUCCESS);
     879             : }
     880             : 
     881             : /*@
     882             :    STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.
     883             : 
     884             :    Not Collective
     885             : 
     886             :    Input Parameter:
     887             : .  st - the spectral transformation context
     888             : 
     889             :    Output Parameters:
     890             : +  m - the number of local rows
     891             : -  n - the number of local columns
     892             : 
     893             :    Level: developer
     894             : 
     895             : .seealso: STMatGetSize()
     896             : @*/
     897         848 : PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)
     898             : {
     899         848 :   PetscFunctionBegin;
     900         848 :   STCheckMatrices(st,1);
     901         848 :   PetscCall(MatGetLocalSize(st->A[0],m,n));
     902         848 :   PetscFunctionReturn(PETSC_SUCCESS);
     903             : }
     904             : 
     905             : /*@C
     906             :    STSetOptionsPrefix - Sets the prefix used for searching for all
     907             :    ST options in the database.
     908             : 
     909             :    Logically Collective
     910             : 
     911             :    Input Parameters:
     912             : +  st     - the spectral transformation context
     913             : -  prefix - the prefix string to prepend to all ST option requests
     914             : 
     915             :    Notes:
     916             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
     917             :    The first character of all runtime options is AUTOMATICALLY the
     918             :    hyphen.
     919             : 
     920             :    Level: advanced
     921             : 
     922             : .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
     923             : @*/
     924         210 : PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)
     925             : {
     926         210 :   PetscFunctionBegin;
     927         210 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     928         210 :   if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp));
     929         210 :   PetscCall(KSPSetOptionsPrefix(st->ksp,prefix));
     930         210 :   PetscCall(KSPAppendOptionsPrefix(st->ksp,"st_"));
     931         210 :   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)st,prefix));
     932         210 :   PetscFunctionReturn(PETSC_SUCCESS);
     933             : }
     934             : 
     935             : /*@C
     936             :    STAppendOptionsPrefix - Appends to the prefix used for searching for all
     937             :    ST options in the database.
     938             : 
     939             :    Logically Collective
     940             : 
     941             :    Input Parameters:
     942             : +  st     - the spectral transformation context
     943             : -  prefix - the prefix string to prepend to all ST option requests
     944             : 
     945             :    Notes:
     946             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
     947             :    The first character of all runtime options is AUTOMATICALLY the
     948             :    hyphen.
     949             : 
     950             :    Level: advanced
     951             : 
     952             : .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
     953             : @*/
     954         174 : PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)
     955             : {
     956         174 :   PetscFunctionBegin;
     957         174 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     958         174 :   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)st,prefix));
     959         174 :   if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp));
     960         174 :   PetscCall(KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix));
     961         174 :   PetscCall(KSPAppendOptionsPrefix(st->ksp,"st_"));
     962         174 :   PetscFunctionReturn(PETSC_SUCCESS);
     963             : }
     964             : 
     965             : /*@C
     966             :    STGetOptionsPrefix - Gets the prefix used for searching for all
     967             :    ST options in the database.
     968             : 
     969             :    Not Collective
     970             : 
     971             :    Input Parameters:
     972             : .  st - the spectral transformation context
     973             : 
     974             :    Output Parameters:
     975             : .  prefix - pointer to the prefix string used, is returned
     976             : 
     977             :    Note:
     978             :    On the Fortran side, the user should pass in a string 'prefix' of
     979             :    sufficient length to hold the prefix.
     980             : 
     981             :    Level: advanced
     982             : 
     983             : .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
     984             : @*/
     985           0 : PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])
     986             : {
     987           0 :   PetscFunctionBegin;
     988           0 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     989           0 :   PetscAssertPointer(prefix,2);
     990           0 :   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)st,prefix));
     991           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     992             : }
     993             : 
     994             : /*@C
     995             :    STView - Prints the ST data structure.
     996             : 
     997             :    Collective
     998             : 
     999             :    Input Parameters:
    1000             : +  st - the ST context
    1001             : -  viewer - optional visualization context
    1002             : 
    1003             :    Note:
    1004             :    The available visualization contexts include
    1005             : +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
    1006             : -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
    1007             :          output where only the first processor opens
    1008             :          the file.  All other processors send their
    1009             :          data to the first processor to print.
    1010             : 
    1011             :    The user can open an alternative visualization contexts with
    1012             :    PetscViewerASCIIOpen() (output to a specified file).
    1013             : 
    1014             :    Level: beginner
    1015             : 
    1016             : .seealso: EPSView()
    1017             : @*/
    1018          11 : PetscErrorCode STView(ST st,PetscViewer viewer)
    1019             : {
    1020          11 :   STType         cstr;
    1021          11 :   char           str[50];
    1022          11 :   PetscBool      isascii,isstring;
    1023             : 
    1024          11 :   PetscFunctionBegin;
    1025          11 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
    1026          11 :   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)st),&viewer));
    1027          11 :   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
    1028          11 :   PetscCheckSameComm(st,1,viewer,2);
    1029             : 
    1030          11 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
    1031          11 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
    1032          11 :   if (isascii) {
    1033          11 :     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer));
    1034          11 :     PetscCall(PetscViewerASCIIPushTab(viewer));
    1035          11 :     PetscTryTypeMethod(st,view,viewer);
    1036          11 :     PetscCall(PetscViewerASCIIPopTab(viewer));
    1037          11 :     PetscCall(SlepcSNPrintfScalar(str,sizeof(str),st->sigma,PETSC_FALSE));
    1038          11 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  shift: %s\n",str));
    1039          11 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  number of matrices: %" PetscInt_FMT "\n",st->nmat));
    1040          11 :     switch (st->matmode) {
    1041             :     case ST_MATMODE_COPY:
    1042             :       break;
    1043           0 :     case ST_MATMODE_INPLACE:
    1044           0 :       PetscCall(PetscViewerASCIIPrintf(viewer,"  shifting the matrix and unshifting at exit\n"));
    1045             :       break;
    1046           1 :     case ST_MATMODE_SHELL:
    1047           1 :       PetscCall(PetscViewerASCIIPrintf(viewer,"  using a shell matrix\n"));
    1048             :       break;
    1049             :     }
    1050          11 :     if (st->nmat>1 && st->matmode != ST_MATMODE_SHELL) PetscCall(PetscViewerASCIIPrintf(viewer,"  nonzero pattern of the matrices: %s\n",MatStructures[st->str]));
    1051          11 :     if (st->Psplit) PetscCall(PetscViewerASCIIPrintf(viewer,"  using split preconditioner matrices with %s\n",MatStructures[st->strp]));
    1052          11 :     if (st->transform && st->nmat>2) PetscCall(PetscViewerASCIIPrintf(viewer,"  computing transformed matrices\n"));
    1053           0 :   } else if (isstring) {
    1054           0 :     PetscCall(STGetType(st,&cstr));
    1055           0 :     PetscCall(PetscViewerStringSPrintf(viewer," %-7.7s",cstr));
    1056           0 :     PetscTryTypeMethod(st,view,viewer);
    1057             :   }
    1058          11 :   if (st->usesksp) {
    1059           6 :     if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp));
    1060           6 :     PetscCall(PetscViewerASCIIPushTab(viewer));
    1061           6 :     PetscCall(KSPView(st->ksp,viewer));
    1062           6 :     PetscCall(PetscViewerASCIIPopTab(viewer));
    1063             :   }
    1064          11 :   PetscFunctionReturn(PETSC_SUCCESS);
    1065             : }
    1066             : 
    1067             : /*@C
    1068             :    STViewFromOptions - View from options
    1069             : 
    1070             :    Collective
    1071             : 
    1072             :    Input Parameters:
    1073             : +  st   - the spectral transformation context
    1074             : .  obj  - optional object
    1075             : -  name - command line option
    1076             : 
    1077             :    Level: intermediate
    1078             : 
    1079             : .seealso: STView(), STCreate()
    1080             : @*/
    1081           0 : PetscErrorCode STViewFromOptions(ST st,PetscObject obj,const char name[])
    1082             : {
    1083           0 :   PetscFunctionBegin;
    1084           0 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
    1085           0 :   PetscCall(PetscObjectViewFromOptions((PetscObject)st,obj,name));
    1086           0 :   PetscFunctionReturn(PETSC_SUCCESS);
    1087             : }
    1088             : 
    1089             : /*@C
    1090             :    STRegister - Adds a method to the spectral transformation package.
    1091             : 
    1092             :    Not Collective
    1093             : 
    1094             :    Input Parameters:
    1095             : +  name - name of a new user-defined transformation
    1096             : -  function - routine to create method context
    1097             : 
    1098             :    Notes:
    1099             :    STRegister() may be called multiple times to add several user-defined
    1100             :    spectral transformations.
    1101             : 
    1102             :    Example Usage:
    1103             : .vb
    1104             :     STRegister("my_transform",MyTransformCreate);
    1105             : .ve
    1106             : 
    1107             :    Then, your spectral transform can be chosen with the procedural interface via
    1108             : $     STSetType(st,"my_transform")
    1109             :    or at runtime via the option
    1110             : $     -st_type my_transform
    1111             : 
    1112             :    Level: advanced
    1113             : 
    1114             : .seealso: STRegisterAll()
    1115             : @*/
    1116        4338 : PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))
    1117             : {
    1118        4338 :   PetscFunctionBegin;
    1119        4338 :   PetscCall(STInitializePackage());
    1120        4338 :   PetscCall(PetscFunctionListAdd(&STList,name,function));
    1121        4338 :   PetscFunctionReturn(PETSC_SUCCESS);
    1122             : }

Generated by: LCOV version 1.14