LCOV - code coverage report
Current view: top level - sys/classes/st/impls/shell - shell.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 101 108 93.5 %
Date: 2024-11-24 00:52:48 Functions: 14 15 93.3 %
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             :    This provides a simple shell interface for programmers to create
      12             :    their own spectral transformations without writing much interface code
      13             : */
      14             : 
      15             : #include <slepc/private/stimpl.h>        /*I "slepcst.h" I*/
      16             : 
      17             : typedef struct {
      18             :   void           *ctx;                       /* user provided context */
      19             :   PetscErrorCode (*apply)(ST,Vec,Vec);
      20             :   PetscErrorCode (*applytrans)(ST,Vec,Vec);
      21             :   PetscErrorCode (*applyhermtrans)(ST,Vec,Vec);
      22             :   PetscErrorCode (*backtransform)(ST,PetscInt n,PetscScalar*,PetscScalar*);
      23             : } ST_SHELL;
      24             : 
      25             : /*@C
      26             :    STShellGetContext - Returns the user-provided context associated with a shell ST
      27             : 
      28             :    Not Collective
      29             : 
      30             :    Input Parameter:
      31             : .  st - spectral transformation context
      32             : 
      33             :    Output Parameter:
      34             : .  ctx - the user provided context
      35             : 
      36             :    Level: advanced
      37             : 
      38             :    Notes:
      39             :    This routine is intended for use within various shell routines
      40             : 
      41             : .seealso: STShellSetContext()
      42             : @*/
      43       22791 : PetscErrorCode STShellGetContext(ST st,void *ctx)
      44             : {
      45       22791 :   PetscBool      flg;
      46             : 
      47       22791 :   PetscFunctionBegin;
      48       22791 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
      49       22791 :   PetscAssertPointer(ctx,2);
      50       22791 :   PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
      51       22791 :   if (!flg) *(void**)ctx = NULL;
      52       22791 :   else      *(void**)ctx = ((ST_SHELL*)st->data)->ctx;
      53       22791 :   PetscFunctionReturn(PETSC_SUCCESS);
      54             : }
      55             : 
      56             : /*@
      57             :    STShellSetContext - Sets the context for a shell ST
      58             : 
      59             :    Logically Collective
      60             : 
      61             :    Input Parameters:
      62             : +  st - the shell ST
      63             : -  ctx - the context
      64             : 
      65             :    Level: advanced
      66             : 
      67             :    Fortran Notes:
      68             :    To use this from Fortran you must write a Fortran interface definition
      69             :    for this function that tells Fortran the Fortran derived data type that
      70             :    you are passing in as the ctx argument.
      71             : 
      72             : .seealso: STShellGetContext()
      73             : @*/
      74          31 : PetscErrorCode STShellSetContext(ST st,void *ctx)
      75             : {
      76          31 :   ST_SHELL       *shell = (ST_SHELL*)st->data;
      77          31 :   PetscBool      flg;
      78             : 
      79          31 :   PetscFunctionBegin;
      80          31 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
      81          31 :   PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
      82          31 :   if (flg) shell->ctx = ctx;
      83          31 :   PetscFunctionReturn(PETSC_SUCCESS);
      84             : }
      85             : 
      86        2842 : static PetscErrorCode STApply_Shell(ST st,Vec x,Vec y)
      87             : {
      88        2842 :   ST_SHELL         *shell = (ST_SHELL*)st->data;
      89        2842 :   PetscObjectState instate,outstate;
      90             : 
      91        2842 :   PetscFunctionBegin;
      92        2842 :   PetscCheck(shell->apply,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No apply() routine provided to Shell ST");
      93        2842 :   PetscCall(VecGetState(y,&instate));
      94        2842 :   PetscCallBack("STSHELL user function apply()",(*shell->apply)(st,x,y));
      95        2842 :   PetscCall(VecGetState(y,&outstate));
      96        2842 :   if (instate == outstate) {
      97             :     /* user forgot to increase the state of the output vector */
      98           0 :     PetscCall(PetscObjectStateIncrease((PetscObject)y));
      99             :   }
     100        2842 :   PetscFunctionReturn(PETSC_SUCCESS);
     101             : }
     102             : 
     103          76 : static PetscErrorCode STApplyTranspose_Shell(ST st,Vec x,Vec y)
     104             : {
     105          76 :   ST_SHELL         *shell = (ST_SHELL*)st->data;
     106          76 :   PetscObjectState instate,outstate;
     107             : 
     108          76 :   PetscFunctionBegin;
     109          76 :   PetscCheck(shell->applytrans,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No applytrans() routine provided to Shell ST");
     110          76 :   PetscCall(VecGetState(y,&instate));
     111          76 :   PetscCallBack("STSHELL user function applytrans()",(*shell->applytrans)(st,x,y));
     112          76 :   PetscCall(VecGetState(y,&outstate));
     113          76 :   if (instate == outstate) {
     114             :     /* user forgot to increase the state of the output vector */
     115           0 :     PetscCall(PetscObjectStateIncrease((PetscObject)y));
     116             :   }
     117          76 :   PetscFunctionReturn(PETSC_SUCCESS);
     118             : }
     119             : 
     120             : #if defined(PETSC_USE_COMPLEX)
     121             : static PetscErrorCode STApplyHermitianTranspose_Shell(ST st,Vec x,Vec y)
     122             : {
     123             :   ST_SHELL         *shell = (ST_SHELL*)st->data;
     124             :   PetscObjectState instate,outstate;
     125             :   Vec              w;
     126             : 
     127             :   PetscFunctionBegin;
     128             :   if (shell->applyhermtrans) {
     129             :     PetscCall(VecGetState(y,&instate));
     130             :     PetscCallBack("STSHELL user function applyhermtrans()",(*shell->applyhermtrans)(st,x,y));
     131             :     PetscCall(VecGetState(y,&outstate));
     132             :     if (instate == outstate) {
     133             :       /* user forgot to increase the state of the output vector */
     134             :       PetscCall(PetscObjectStateIncrease((PetscObject)y));
     135             :     }
     136             :   } else {
     137             :     PetscCall(VecDuplicate(x,&w));
     138             :     PetscCall(VecCopy(x,w));
     139             :     PetscCall(VecConjugate(w));
     140             :     PetscCall(STApplyTranspose_Shell(st,w,y));
     141             :     PetscCall(VecDestroy(&w));
     142             :     PetscCall(VecConjugate(y));
     143             :   }
     144             :   PetscFunctionReturn(PETSC_SUCCESS);
     145             : }
     146             : #endif
     147             : 
     148       21466 : static PetscErrorCode STBackTransform_Shell(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
     149             : {
     150       21466 :   ST_SHELL       *shell = (ST_SHELL*)st->data;
     151             : 
     152       21466 :   PetscFunctionBegin;
     153       21466 :   if (shell->backtransform) PetscCallBack("STSHELL user function backtransform()",(*shell->backtransform)(st,n,eigr,eigi));
     154       21466 :   PetscFunctionReturn(PETSC_SUCCESS);
     155             : }
     156             : 
     157             : /*
     158             :    STIsInjective_Shell - Check if the user has provided the backtransform operation.
     159             : */
     160          31 : PetscErrorCode STIsInjective_Shell(ST st,PetscBool* is)
     161             : {
     162          31 :   ST_SHELL *shell = (ST_SHELL*)st->data;
     163             : 
     164          31 :   PetscFunctionBegin;
     165          31 :   *is = shell->backtransform? PETSC_TRUE: PETSC_FALSE;
     166          31 :   PetscFunctionReturn(PETSC_SUCCESS);
     167             : }
     168             : 
     169          35 : static PetscErrorCode STDestroy_Shell(ST st)
     170             : {
     171          35 :   PetscFunctionBegin;
     172          35 :   PetscCall(PetscFree(st->data));
     173          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",NULL));
     174          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",NULL));
     175          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",NULL));
     176          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",NULL));
     177          35 :   PetscFunctionReturn(PETSC_SUCCESS);
     178             : }
     179             : 
     180          31 : static PetscErrorCode STShellSetApply_Shell(ST st,PetscErrorCode (*apply)(ST,Vec,Vec))
     181             : {
     182          31 :   ST_SHELL *shell = (ST_SHELL*)st->data;
     183             : 
     184          31 :   PetscFunctionBegin;
     185          31 :   shell->apply = apply;
     186          31 :   PetscFunctionReturn(PETSC_SUCCESS);
     187             : }
     188             : 
     189             : /*@C
     190             :    STShellSetApply - Sets routine to use as the application of the
     191             :    operator to a vector in the user-defined spectral transformation.
     192             : 
     193             :    Logically Collective
     194             : 
     195             :    Input Parameters:
     196             : +  st    - the spectral transformation context
     197             : -  apply - the application-provided transformation routine
     198             : 
     199             :    Calling sequence of apply:
     200             : $  PetscErrorCode apply(ST st,Vec xin,Vec xout)
     201             : +  st   - the spectral transformation context
     202             : .  xin  - input vector
     203             : -  xout - output vector
     204             : 
     205             :    Level: advanced
     206             : 
     207             : .seealso: STShellSetBackTransform(), STShellSetApplyTranspose(), STShellSetApplyHermitianTranspose()
     208             : @*/
     209          31 : PetscErrorCode STShellSetApply(ST st,PetscErrorCode (*apply)(ST st,Vec xin,Vec xout))
     210             : {
     211          31 :   PetscFunctionBegin;
     212          31 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     213          31 :   PetscTryMethod(st,"STShellSetApply_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,apply));
     214          31 :   PetscFunctionReturn(PETSC_SUCCESS);
     215             : }
     216             : 
     217           6 : static PetscErrorCode STShellSetApplyTranspose_Shell(ST st,PetscErrorCode (*applytrans)(ST,Vec,Vec))
     218             : {
     219           6 :   ST_SHELL *shell = (ST_SHELL*)st->data;
     220             : 
     221           6 :   PetscFunctionBegin;
     222           6 :   shell->applytrans = applytrans;
     223           6 :   PetscFunctionReturn(PETSC_SUCCESS);
     224             : }
     225             : 
     226             : /*@C
     227             :    STShellSetApplyTranspose - Sets routine to use as the application of the
     228             :    transposed operator to a vector in the user-defined spectral transformation.
     229             : 
     230             :    Logically Collective
     231             : 
     232             :    Input Parameters:
     233             : +  st    - the spectral transformation context
     234             : -  applytrans - the application-provided transformation routine
     235             : 
     236             :    Calling sequence of applytrans:
     237             : $  PetscErrorCode applytrans(ST st,Vec xin,Vec xout)
     238             : +  st   - the spectral transformation context
     239             : .  xin  - input vector
     240             : -  xout - output vector
     241             : 
     242             :    Level: advanced
     243             : 
     244             : .seealso: STShellSetApply(), STShellSetBackTransform()
     245             : @*/
     246           6 : PetscErrorCode STShellSetApplyTranspose(ST st,PetscErrorCode (*applytrans)(ST st,Vec xin,Vec xout))
     247             : {
     248           6 :   PetscFunctionBegin;
     249           6 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     250           6 :   PetscTryMethod(st,"STShellSetApplyTranspose_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,applytrans));
     251           6 :   PetscFunctionReturn(PETSC_SUCCESS);
     252             : }
     253             : 
     254             : #if defined(PETSC_USE_COMPLEX)
     255             : static PetscErrorCode STShellSetApplyHermitianTranspose_Shell(ST st,PetscErrorCode (*applyhermtrans)(ST,Vec,Vec))
     256             : {
     257             :   ST_SHELL *shell = (ST_SHELL*)st->data;
     258             : 
     259             :   PetscFunctionBegin;
     260             :   shell->applyhermtrans = applyhermtrans;
     261             :   PetscFunctionReturn(PETSC_SUCCESS);
     262             : }
     263             : #endif
     264             : 
     265             : /*@C
     266             :    STShellSetApplyHermitianTranspose - Sets routine to use as the application of the
     267             :    conjugate-transposed operator to a vector in the user-defined spectral transformation.
     268             : 
     269             :    Logically Collective
     270             : 
     271             :    Input Parameters:
     272             : +  st    - the spectral transformation context
     273             : -  applyhermtrans - the application-provided transformation routine
     274             : 
     275             :    Calling sequence of applyhermtrans:
     276             : $  PetscErrorCode applyhermtrans(ST st,Vec xin,Vec xout)
     277             : +  st   - the spectral transformation context
     278             : .  xin  - input vector
     279             : -  xout - output vector
     280             : 
     281             :    Note:
     282             :    If configured with real scalars, this function has the same effect as STShellSetApplyTranspose(),
     283             :    so no need to call both.
     284             : 
     285             :    Level: advanced
     286             : 
     287             : .seealso: STShellSetApply(), STShellSetApplyTranspose(), STShellSetBackTransform()
     288             : @*/
     289           0 : PetscErrorCode STShellSetApplyHermitianTranspose(ST st,PetscErrorCode (*applyhermtrans)(ST st,Vec xin,Vec xout))
     290             : {
     291           0 :   PetscFunctionBegin;
     292           0 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     293           0 :   PetscTryMethod(st,"STShellSetApplyHermitianTranspose_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,applyhermtrans));
     294           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     295             : }
     296             : 
     297          30 : static PetscErrorCode STShellSetBackTransform_Shell(ST st,PetscErrorCode (*backtr)(ST,PetscInt,PetscScalar*,PetscScalar*))
     298             : {
     299          30 :   ST_SHELL *shell = (ST_SHELL*)st->data;
     300             : 
     301          30 :   PetscFunctionBegin;
     302          30 :   shell->backtransform = backtr;
     303          30 :   PetscFunctionReturn(PETSC_SUCCESS);
     304             : }
     305             : 
     306             : /*@C
     307             :    STShellSetBackTransform - Sets the routine to be called after the
     308             :    eigensolution process has finished in order to transform back the
     309             :    computed eigenvalues.
     310             : 
     311             :    Logically Collective
     312             : 
     313             :    Input Parameters:
     314             : +  st     - the spectral transformation context
     315             : -  backtr - the application-provided backtransform routine
     316             : 
     317             :    Calling sequence of backtr:
     318             : $  PetscErrorCode backtr(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
     319             : +  st   - the spectral transformation context
     320             : .  n    - number of eigenvalues to be backtransformed
     321             : .  eigr - pointer ot the real parts of the eigenvalues to transform back
     322             : -  eigi - pointer ot the imaginary parts
     323             : 
     324             :    Level: advanced
     325             : 
     326             : .seealso: STShellSetApply(), STShellSetApplyTranspose()
     327             : @*/
     328          30 : PetscErrorCode STShellSetBackTransform(ST st,PetscErrorCode (*backtr)(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi))
     329             : {
     330          30 :   PetscFunctionBegin;
     331          30 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     332          30 :   PetscTryMethod(st,"STShellSetBackTransform_C",(ST,PetscErrorCode (*)(ST,PetscInt,PetscScalar*,PetscScalar*)),(st,backtr));
     333          30 :   PetscFunctionReturn(PETSC_SUCCESS);
     334             : }
     335             : 
     336             : /*MC
     337             :    STSHELL - User-defined spectral transformation via callback functions
     338             :    for the application of the operator to a vector and (optionally) the
     339             :    backtransform operation.
     340             : 
     341             :    Level: advanced
     342             : 
     343             :    Usage:
     344             : $             extern PetscErrorCode (*apply)(void*,Vec,Vec);
     345             : $             extern PetscErrorCode (*applytrans)(void*,Vec,Vec);
     346             : $             extern PetscErrorCode (*applyht)(void*,Vec,Vec);
     347             : $             extern PetscErrorCode (*backtr)(void*,PetscScalar*,PetscScalar*);
     348             : $
     349             : $             STCreate(comm,&st);
     350             : $             STSetType(st,STSHELL);
     351             : $             STShellSetContext(st,ctx);
     352             : $             STShellSetApply(st,apply);
     353             : $             STShellSetApplyTranspose(st,applytrans);        (optional)
     354             : $             STShellSetApplyHermitianTranspose(st,applyht);  (optional, only in complex scalars)
     355             : $             STShellSetBackTransform(st,backtr);             (optional)
     356             : 
     357             : M*/
     358             : 
     359          35 : SLEPC_EXTERN PetscErrorCode STCreate_Shell(ST st)
     360             : {
     361          35 :   ST_SHELL       *ctx;
     362             : 
     363          35 :   PetscFunctionBegin;
     364          35 :   PetscCall(PetscNew(&ctx));
     365          35 :   st->data = (void*)ctx;
     366             : 
     367          35 :   st->usesksp = PETSC_FALSE;
     368             : 
     369          35 :   st->ops->apply           = STApply_Shell;
     370          35 :   st->ops->applytrans      = STApplyTranspose_Shell;
     371             : #if defined(PETSC_USE_COMPLEX)
     372             :   st->ops->applyhermtrans  = STApplyHermitianTranspose_Shell;
     373             : #else
     374          35 :   st->ops->applyhermtrans  = STApplyTranspose_Shell;
     375             : #endif
     376          35 :   st->ops->backtransform   = STBackTransform_Shell;
     377          35 :   st->ops->destroy         = STDestroy_Shell;
     378             : 
     379          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",STShellSetApply_Shell));
     380          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",STShellSetApplyTranspose_Shell));
     381             : #if defined(PETSC_USE_COMPLEX)
     382             :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",STShellSetApplyHermitianTranspose_Shell));
     383             : #else
     384          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",STShellSetApplyTranspose_Shell));
     385             : #endif
     386          35 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",STShellSetBackTransform_Shell));
     387          35 :   PetscFunctionReturn(PETSC_SUCCESS);
     388             : }

Generated by: LCOV version 1.14