Actual source code: shift.c

slepc-3.15.0 2021-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Shift spectral transformation, applies (A + sigma I) as operator, or
 12:    inv(B)(A + sigma B) for generalized problems
 13: */

 15: #include <slepc/private/stimpl.h>

 17: PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 18: {
 19:   PetscInt j;

 22:   for (j=0;j<n;j++) {
 23:     eigr[j] += st->sigma;
 24:   }
 25:   return(0);
 26: }

 28: PetscErrorCode STPostSolve_Shift(ST st)
 29: {

 33:   if (st->matmode == ST_MATMODE_INPLACE) {
 34:     if (st->nmat>1) {
 35:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
 36:     } else {
 37:       MatShift(st->A[0],st->sigma);
 38:     }
 39:     st->Astate[0] = ((PetscObject)st->A[0])->state;
 40:     st->state   = ST_STATE_INITIAL;
 41:     st->opready = PETSC_FALSE;
 42:   }
 43:   return(0);
 44: }

 46: /*
 47:    Operator (shift):
 48:                Op               P         M
 49:    if nmat=1:  A-sI             NULL      A-sI
 50:    if nmat=2:  B^-1 (A-sB)      B         A-sB
 51: */
 52: PetscErrorCode STComputeOperator_Shift(ST st)
 53: {

 57:   st->usesksp = (st->nmat>1)? PETSC_TRUE: PETSC_FALSE;
 58:   PetscObjectReference((PetscObject)st->A[1]);
 59:   MatDestroy(&st->T[1]);
 60:   st->T[1] = st->A[1];
 61:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[0]);
 62:   if (st->nmat>1) { PetscObjectReference((PetscObject)st->T[1]); }
 63:   MatDestroy(&st->P);
 64:   st->P = (st->nmat>1)? st->T[1]: NULL;
 65:   st->M = st->T[0];
 66:   return(0);
 67: }

 69: PetscErrorCode STSetUp_Shift(ST st)
 70: {
 72:   PetscInt       k,nc,nmat=st->nmat;
 73:   PetscScalar    *coeffs=NULL;

 76:   if (nmat>1) {
 77:     STSetWorkVecs(st,1);
 78:   }
 79:   if (nmat>2) {  /* set-up matrices for polynomial eigenproblems */
 80:     if (st->transform) {
 81:       nc = (nmat*(nmat+1))/2;
 82:       PetscMalloc1(nc,&coeffs);
 83:       /* Compute coeffs */
 84:       STCoeffs_Monomial(st,coeffs);
 85:       /* T[n] = A_n */
 86:       k = nmat-1;
 87:       PetscObjectReference((PetscObject)st->A[k]);
 88:       MatDestroy(&st->T[k]);
 89:       st->T[k] = st->A[k];
 90:       for (k=0;k<nmat-1;k++) {
 91:         STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[k]);
 92:       }
 93:       PetscFree(coeffs);
 94:       PetscObjectReference((PetscObject)st->T[nmat-1]);
 95:       MatDestroy(&st->P);
 96:       st->P = st->T[nmat-1];
 97:       STKSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
 98:     } else {
 99:       for (k=0;k<nmat;k++) {
100:         PetscObjectReference((PetscObject)st->A[k]);
101:         MatDestroy(&st->T[k]);
102:         st->T[k] = st->A[k];
103:       }
104:     }
105:   }
106:   if (st->P) {
107:     KSPSetUp(st->ksp);
108:   }
109:   return(0);
110: }

112: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
113: {
115:   PetscInt       k,nc,nmat=PetscMax(st->nmat,2);
116:   PetscScalar    *coeffs=NULL;

119:   if (st->transform) {
120:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
121:       nc = (nmat*(nmat+1))/2;
122:       PetscMalloc1(nc,&coeffs);
123:       /* Compute coeffs */
124:       STCoeffs_Monomial(st,coeffs);
125:     }
126:     for (k=0;k<nmat-1;k++) {
127:       STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PETSC_FALSE,&st->T[k]);
128:     }
129:     if (st->matmode == ST_MATMODE_COPY && nmat>2) {
130:       PetscFree(coeffs);
131:     }
132:     if (st->nmat<=2) st->M = st->T[0];
133:   }
134:   return(0);
135: }

137: SLEPC_EXTERN PetscErrorCode STCreate_Shift(ST st)
138: {
140:   st->usesksp = PETSC_TRUE;

142:   st->ops->apply           = STApply_Generic;
143:   st->ops->applytrans      = STApplyTranspose_Generic;
144:   st->ops->backtransform   = STBackTransform_Shift;
145:   st->ops->setshift        = STSetShift_Shift;
146:   st->ops->getbilinearform = STGetBilinearForm_Default;
147:   st->ops->setup           = STSetUp_Shift;
148:   st->ops->computeoperator = STComputeOperator_Shift;
149:   st->ops->postsolve       = STPostSolve_Shift;
150:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
151:   return(0);
152: }