Actual source code: shift.c
slepc-3.13.2 2020-05-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: if (!st->ksp) { STGetKSP(st,&st->ksp); }
98: STCheckFactorPackage(st);
99: KSPSetOperators(st->ksp,st->P,st->P);
100: } else {
101: for (k=0;k<nmat;k++) {
102: PetscObjectReference((PetscObject)st->A[k]);
103: MatDestroy(&st->T[k]);
104: st->T[k] = st->A[k];
105: }
106: }
107: }
108: if (st->P) {
109: KSPSetUp(st->ksp);
110: }
111: return(0);
112: }
114: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
115: {
117: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
118: PetscScalar *coeffs=NULL;
121: if (st->transform) {
122: if (st->matmode == ST_MATMODE_COPY && nmat>2) {
123: nc = (nmat*(nmat+1))/2;
124: PetscMalloc1(nc,&coeffs);
125: /* Compute coeffs */
126: STCoeffs_Monomial(st,coeffs);
127: }
128: for (k=0;k<nmat-1;k++) {
129: 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]);
130: }
131: if (st->matmode == ST_MATMODE_COPY && nmat>2) {
132: PetscFree(coeffs);
133: }
134: if (st->nmat<=2) st->M = st->T[0];
135: }
136: return(0);
137: }
139: SLEPC_EXTERN PetscErrorCode STCreate_Shift(ST st)
140: {
142: st->usesksp = PETSC_TRUE;
144: st->ops->apply = STApply_Generic;
145: st->ops->applytrans = STApplyTranspose_Generic;
146: st->ops->backtransform = STBackTransform_Shift;
147: st->ops->setshift = STSetShift_Shift;
148: st->ops->getbilinearform = STGetBilinearForm_Default;
149: st->ops->setup = STSetUp_Shift;
150: st->ops->computeoperator = STComputeOperator_Shift;
151: st->ops->postsolve = STPostSolve_Shift;
152: st->ops->setdefaultksp = STSetDefaultKSP_Default;
153: return(0);
154: }