LCOV - code coverage report
Current view: top level - sys/classes/st/impls/cayley - cayley.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 163 176 92.6 %
Date: 2024-03-29 00:34:52 Functions: 15 16 93.8 %
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             :    Implements the Cayley spectral transform
      12             : */
      13             : 
      14             : #include <slepc/private/stimpl.h>          /*I "slepcst.h" I*/
      15             : 
      16             : typedef struct {
      17             :   PetscScalar nu;
      18             :   PetscBool   nu_set;
      19             : } ST_CAYLEY;
      20             : 
      21          16 : static PetscErrorCode MatMult_Cayley(Mat B,Vec x,Vec y)
      22             : {
      23          16 :   ST             st;
      24          16 :   ST_CAYLEY      *ctx;
      25          16 :   PetscScalar    nu;
      26             : 
      27          16 :   PetscFunctionBegin;
      28          16 :   PetscCall(MatShellGetContext(B,&st));
      29          16 :   ctx = (ST_CAYLEY*)st->data;
      30          16 :   nu = ctx->nu;
      31             : 
      32          16 :   if (st->matmode == ST_MATMODE_INPLACE) { nu = nu + st->sigma; }
      33             : 
      34          16 :   if (st->nmat>1) {
      35             :     /* generalized eigenproblem: y = (A + tB)x */
      36           4 :     PetscCall(MatMult(st->A[0],x,y));
      37           4 :     PetscCall(MatMult(st->A[1],x,st->work[1]));
      38           4 :     PetscCall(VecAXPY(y,nu,st->work[1]));
      39             :   } else {
      40             :     /* standard eigenproblem: y = (A + tI)x */
      41          12 :     PetscCall(MatMult(st->A[0],x,y));
      42          12 :     PetscCall(VecAXPY(y,nu,x));
      43             :   }
      44          16 :   PetscFunctionReturn(PETSC_SUCCESS);
      45             : }
      46             : 
      47           6 : static PetscErrorCode MatMultTranspose_Cayley(Mat B,Vec x,Vec y)
      48             : {
      49           6 :   ST             st;
      50           6 :   ST_CAYLEY      *ctx;
      51           6 :   PetscScalar    nu;
      52             : 
      53           6 :   PetscFunctionBegin;
      54           6 :   PetscCall(MatShellGetContext(B,&st));
      55           6 :   ctx = (ST_CAYLEY*)st->data;
      56           6 :   nu = ctx->nu;
      57             : 
      58           6 :   if (st->matmode == ST_MATMODE_INPLACE) { nu = nu + st->sigma; }
      59           6 :   nu = PetscConj(nu);
      60             : 
      61           6 :   if (st->nmat>1) {
      62             :     /* generalized eigenproblem: y = (A + tB)x */
      63           0 :     PetscCall(MatMultTranspose(st->A[0],x,y));
      64           0 :     PetscCall(MatMultTranspose(st->A[1],x,st->work[1]));
      65           0 :     PetscCall(VecAXPY(y,nu,st->work[1]));
      66             :   } else {
      67             :     /* standard eigenproblem: y = (A + tI)x */
      68           6 :     PetscCall(MatMultTranspose(st->A[0],x,y));
      69           6 :     PetscCall(VecAXPY(y,nu,x));
      70             :   }
      71           6 :   PetscFunctionReturn(PETSC_SUCCESS);
      72             : }
      73             : 
      74          10 : static PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
      75             : {
      76          10 :   PetscFunctionBegin;
      77          10 :   PetscCall(STSetUp(st));
      78          10 :   *B = st->T[0];
      79          10 :   PetscCall(PetscObjectReference((PetscObject)*B));
      80          10 :   PetscFunctionReturn(PETSC_SUCCESS);
      81             : }
      82             : 
      83         745 : static PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
      84             : {
      85         745 :   ST_CAYLEY   *ctx = (ST_CAYLEY*)st->data;
      86         745 :   PetscInt    j;
      87             : #if !defined(PETSC_USE_COMPLEX)
      88             :   PetscScalar t,i,r;
      89             : #endif
      90             : 
      91         745 :   PetscFunctionBegin;
      92             : #if !defined(PETSC_USE_COMPLEX)
      93             :   for (j=0;j<n;j++) {
      94             :     if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
      95             :     else {
      96             :       r = eigr[j];
      97             :       i = eigi[j];
      98             :       r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
      99             :       i = - st->sigma * i - ctx->nu * i;
     100             :       t = i * i + r * (r - 2.0) + 1.0;
     101             :       eigr[j] = r / t;
     102             :       eigi[j] = i / t;
     103             :     }
     104             :   }
     105             : #else
     106        2233 :   for (j=0;j<n;j++) {
     107        1488 :     eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
     108             :   }
     109             : #endif
     110         745 :   PetscFunctionReturn(PETSC_SUCCESS);
     111             : }
     112             : 
     113           8 : static PetscErrorCode STPostSolve_Cayley(ST st)
     114             : {
     115           8 :   PetscFunctionBegin;
     116           8 :   if (st->matmode == ST_MATMODE_INPLACE) {
     117           2 :     if (st->nmat>1) PetscCall(MatAXPY(st->A[0],st->sigma,st->A[1],st->str));
     118           2 :     else PetscCall(MatShift(st->A[0],st->sigma));
     119           2 :     st->Astate[0] = ((PetscObject)st->A[0])->state;
     120           2 :     st->state   = ST_STATE_INITIAL;
     121           2 :     st->opready = PETSC_FALSE;
     122             :   }
     123           8 :   PetscFunctionReturn(PETSC_SUCCESS);
     124             : }
     125             : 
     126             : /*
     127             :    Operator (cayley):
     128             :                Op                  P         M
     129             :    if nmat=1:  (A-sI)^-1 (A+tI)    A-sI      A+tI
     130             :    if nmat=2:  (A-sB)^-1 (A+tB)    A-sB      A+tI
     131             : */
     132          19 : static PetscErrorCode STComputeOperator_Cayley(ST st)
     133             : {
     134          19 :   PetscInt       n,m;
     135          19 :   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
     136             : 
     137          19 :   PetscFunctionBegin;
     138             :   /* if the user did not set the shift, use the target value */
     139          19 :   if (!st->sigma_set) st->sigma = st->defsigma;
     140             : 
     141          19 :   if (!ctx->nu_set) ctx->nu = st->sigma;
     142          19 :   PetscCheck(ctx->nu!=0.0 || st->sigma!=0.0,PetscObjectComm((PetscObject)st),PETSC_ERR_USER_INPUT,"Values of shift and antishift cannot be zero simultaneously");
     143          19 :   PetscCheck(ctx->nu!=-st->sigma,PetscObjectComm((PetscObject)st),PETSC_ERR_USER_INPUT,"It is not allowed to set the antishift equal to minus the shift (the target)");
     144             : 
     145             :   /* T[0] = A+nu*B */
     146          19 :   if (st->matmode==ST_MATMODE_INPLACE) {
     147           5 :     PetscCall(MatGetLocalSize(st->A[0],&n,&m));
     148           5 :     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)st),n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,&st->T[0]));
     149           5 :     PetscCall(MatShellSetOperation(st->T[0],MATOP_MULT,(void(*)(void))MatMult_Cayley));
     150           5 :     PetscCall(MatShellSetOperation(st->T[0],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Cayley));
     151          14 :   } else PetscCall(STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[0]));
     152          19 :   st->M = st->T[0];
     153             : 
     154             :   /* T[1] = A-sigma*B */
     155          19 :   PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]));
     156          19 :   PetscCall(PetscObjectReference((PetscObject)st->T[1]));
     157          19 :   PetscCall(MatDestroy(&st->P));
     158          19 :   st->P = st->T[1];
     159          19 :   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
     160           2 :     PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
     161             :   }
     162          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     163             : }
     164             : 
     165          19 : static PetscErrorCode STSetUp_Cayley(ST st)
     166             : {
     167          19 :   PetscFunctionBegin;
     168          19 :   PetscCheck(st->nmat<=2,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cayley transform cannot be used in polynomial eigenproblems");
     169          19 :   PetscCall(STSetWorkVecs(st,2));
     170          19 :   PetscCall(KSPSetUp(st->ksp));
     171          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     172             : }
     173             : 
     174          11 : static PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
     175             : {
     176          11 :   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
     177             : 
     178          11 :   PetscFunctionBegin;
     179          11 :   PetscCheck(newshift!=0.0 || (ctx->nu_set && ctx->nu!=0.0),PetscObjectComm((PetscObject)st),PETSC_ERR_USER_INPUT,"Values of shift and antishift cannot be zero simultaneously");
     180          11 :   PetscCheck(ctx->nu!=-newshift,PetscObjectComm((PetscObject)st),PETSC_ERR_USER_INPUT,"It is not allowed to set the shift equal to minus the antishift");
     181             : 
     182          11 :   if (!ctx->nu_set) {
     183           9 :     if (st->matmode!=ST_MATMODE_INPLACE) PetscCall(STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]));
     184           9 :     ctx->nu = newshift;
     185             :   }
     186          11 :   PetscCall(STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[1]));
     187          11 :   if (st->P!=st->T[1]) {
     188           2 :     PetscCall(PetscObjectReference((PetscObject)st->T[1]));
     189           2 :     PetscCall(MatDestroy(&st->P));
     190           2 :     st->P = st->T[1];
     191             :   }
     192          11 :   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
     193           2 :     PetscCall(STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_TRUE,&st->Pmat));
     194             :   }
     195          11 :   PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
     196          11 :   PetscCall(KSPSetUp(st->ksp));
     197          11 :   PetscFunctionReturn(PETSC_SUCCESS);
     198             : }
     199             : 
     200          10 : static PetscErrorCode STSetFromOptions_Cayley(ST st,PetscOptionItems *PetscOptionsObject)
     201             : {
     202          10 :   PetscScalar    nu;
     203          10 :   PetscBool      flg;
     204          10 :   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
     205             : 
     206          10 :   PetscFunctionBegin;
     207          10 :   PetscOptionsHeadBegin(PetscOptionsObject,"ST Cayley Options");
     208             : 
     209          10 :     PetscCall(PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg));
     210          10 :     if (flg) PetscCall(STCayleySetAntishift(st,nu));
     211             : 
     212          10 :   PetscOptionsHeadEnd();
     213          10 :   PetscFunctionReturn(PETSC_SUCCESS);
     214             : }
     215             : 
     216          18 : static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
     217             : {
     218          18 :   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
     219             : 
     220          18 :   PetscFunctionBegin;
     221          18 :   if (ctx->nu != newshift) {
     222          18 :     STCheckNotSeized(st,1);
     223          18 :     if (st->state && st->matmode!=ST_MATMODE_INPLACE) PetscCall(STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]));
     224          18 :     ctx->nu = newshift;
     225             :   }
     226          18 :   ctx->nu_set = PETSC_TRUE;
     227          18 :   PetscFunctionReturn(PETSC_SUCCESS);
     228             : }
     229             : 
     230             : /*@
     231             :    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
     232             :    spectral transformation.
     233             : 
     234             :    Logically Collective
     235             : 
     236             :    Input Parameters:
     237             : +  st  - the spectral transformation context
     238             : -  nu  - the anti-shift
     239             : 
     240             :    Options Database Key:
     241             : .  -st_cayley_antishift - Sets the value of the anti-shift
     242             : 
     243             :    Level: intermediate
     244             : 
     245             :    Note:
     246             :    In the generalized Cayley transform, the operator can be expressed as
     247             :    OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
     248             :    Use STSetShift() for setting sigma. The value nu=-sigma is not allowed.
     249             : 
     250             : .seealso: STSetShift(), STCayleyGetAntishift()
     251             : @*/
     252          22 : PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
     253             : {
     254          22 :   PetscFunctionBegin;
     255          22 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     256          88 :   PetscValidLogicalCollectiveScalar(st,nu,2);
     257          22 :   PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
     258          22 :   PetscFunctionReturn(PETSC_SUCCESS);
     259             : }
     260             : 
     261          27 : static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
     262             : {
     263          27 :   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
     264             : 
     265          27 :   PetscFunctionBegin;
     266          27 :   *nu = ctx->nu;
     267          27 :   PetscFunctionReturn(PETSC_SUCCESS);
     268             : }
     269             : 
     270             : /*@
     271             :    STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
     272             :    spectral transformation.
     273             : 
     274             :    Not Collective
     275             : 
     276             :    Input Parameter:
     277             : .  st  - the spectral transformation context
     278             : 
     279             :    Output Parameter:
     280             : .  nu  - the anti-shift
     281             : 
     282             :    Level: intermediate
     283             : 
     284             : .seealso: STGetShift(), STCayleySetAntishift()
     285             : @*/
     286          27 : PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
     287             : {
     288          27 :   PetscFunctionBegin;
     289          27 :   PetscValidHeaderSpecific(st,ST_CLASSID,1);
     290          27 :   PetscAssertPointer(nu,2);
     291          27 :   PetscUseMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
     292          27 :   PetscFunctionReturn(PETSC_SUCCESS);
     293             : }
     294             : 
     295           0 : static PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
     296             : {
     297           0 :   char           str[50];
     298           0 :   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
     299           0 :   PetscBool      isascii;
     300             : 
     301           0 :   PetscFunctionBegin;
     302           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     303           0 :   if (isascii) {
     304           0 :     PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->nu,PETSC_FALSE));
     305           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  antishift: %s\n",str));
     306             :   }
     307           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     308             : }
     309             : 
     310          19 : static PetscErrorCode STDestroy_Cayley(ST st)
     311             : {
     312          19 :   PetscFunctionBegin;
     313          19 :   PetscCall(PetscFree(st->data));
     314          19 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL));
     315          19 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL));
     316          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     317             : }
     318             : 
     319          19 : SLEPC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
     320             : {
     321          19 :   ST_CAYLEY      *ctx;
     322             : 
     323          19 :   PetscFunctionBegin;
     324          19 :   PetscCall(PetscNew(&ctx));
     325          19 :   st->data = (void*)ctx;
     326             : 
     327          19 :   st->usesksp = PETSC_TRUE;
     328             : 
     329          19 :   st->ops->apply           = STApply_Generic;
     330          19 :   st->ops->applytrans      = STApplyTranspose_Generic;
     331          19 :   st->ops->applyhermtrans  = STApplyHermitianTranspose_Generic;
     332          19 :   st->ops->backtransform   = STBackTransform_Cayley;
     333          19 :   st->ops->setshift        = STSetShift_Cayley;
     334          19 :   st->ops->getbilinearform = STGetBilinearForm_Cayley;
     335          19 :   st->ops->setup           = STSetUp_Cayley;
     336          19 :   st->ops->computeoperator = STComputeOperator_Cayley;
     337          19 :   st->ops->setfromoptions  = STSetFromOptions_Cayley;
     338          19 :   st->ops->postsolve       = STPostSolve_Cayley;
     339          19 :   st->ops->destroy         = STDestroy_Cayley;
     340          19 :   st->ops->view            = STView_Cayley;
     341          19 :   st->ops->checknullspace  = STCheckNullSpace_Default;
     342          19 :   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
     343             : 
     344          19 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley));
     345          19 :   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley));
     346          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     347             : }

Generated by: LCOV version 1.14