Actual source code: shift.c
slepc-3.21.1 2024-04-26
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: static PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
18: {
19: PetscInt j;
21: PetscFunctionBegin;
22: for (j=0;j<n;j++) {
23: eigr[j] += st->sigma;
24: }
25: PetscFunctionReturn(PETSC_SUCCESS);
26: }
28: static PetscErrorCode STPostSolve_Shift(ST st)
29: {
30: PetscFunctionBegin;
31: if (st->matmode == ST_MATMODE_INPLACE) {
32: if (st->nmat>1) PetscCall(MatAXPY(st->A[0],st->sigma,st->A[1],st->str));
33: else PetscCall(MatShift(st->A[0],st->sigma));
34: st->Astate[0] = ((PetscObject)st->A[0])->state;
35: st->state = ST_STATE_INITIAL;
36: st->opready = PETSC_FALSE;
37: }
38: PetscFunctionReturn(PETSC_SUCCESS);
39: }
41: /*
42: Operator (shift):
43: Op P M
44: if nmat=1: A-sI NULL A-sI
45: if nmat=2: B^-1 (A-sB) B A-sB
46: */
47: static PetscErrorCode STComputeOperator_Shift(ST st)
48: {
49: PetscFunctionBegin;
50: st->usesksp = (st->nmat>1)? PETSC_TRUE: PETSC_FALSE;
51: PetscCall(PetscObjectReference((PetscObject)st->A[1]));
52: PetscCall(MatDestroy(&st->T[1]));
53: st->T[1] = st->A[1];
54: PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[0]));
55: if (st->nmat>1) PetscCall(PetscObjectReference((PetscObject)st->T[1]));
56: PetscCall(MatDestroy(&st->P));
57: st->P = (st->nmat>1)? st->T[1]: NULL;
58: st->M = st->T[0];
59: if (st->nmat>1 && st->Psplit) { /* build custom preconditioner from the split matrices */
60: PetscCall(MatDestroy(&st->Pmat));
61: PetscCall(PetscObjectReference((PetscObject)st->Psplit[1]));
62: st->Pmat = st->Psplit[1];
63: }
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: static PetscErrorCode STSetUp_Shift(ST st)
68: {
69: PetscInt k,nc,nmat=st->nmat;
70: PetscScalar *coeffs=NULL;
72: PetscFunctionBegin;
73: if (nmat>1) PetscCall(STSetWorkVecs(st,1));
74: if (nmat>2) { /* set-up matrices for polynomial eigenproblems */
75: if (st->transform) {
76: nc = (nmat*(nmat+1))/2;
77: PetscCall(PetscMalloc1(nc,&coeffs));
78: /* Compute coeffs */
79: PetscCall(STCoeffs_Monomial(st,coeffs));
80: /* T[n] = A_n */
81: k = nmat-1;
82: PetscCall(PetscObjectReference((PetscObject)st->A[k]));
83: PetscCall(MatDestroy(&st->T[k]));
84: st->T[k] = st->A[k];
85: for (k=0;k<nmat-1;k++) PetscCall(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),PETSC_FALSE,&st->T[k]));
86: PetscCall(PetscFree(coeffs));
87: PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
88: PetscCall(MatDestroy(&st->P));
89: st->P = st->T[nmat-1];
90: if (st->Psplit) { /* build custom preconditioner from the split matrices */
91: PetscCall(STMatMAXPY_Private(st,st->sigma,0.0,nmat-1,coeffs?coeffs:NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
92: }
93: PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
94: } else {
95: for (k=0;k<nmat;k++) {
96: PetscCall(PetscObjectReference((PetscObject)st->A[k]));
97: PetscCall(MatDestroy(&st->T[k]));
98: st->T[k] = st->A[k];
99: }
100: }
101: }
102: if (st->P) PetscCall(KSPSetUp(st->ksp));
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: static PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
107: {
108: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
109: PetscScalar *coeffs=NULL;
111: PetscFunctionBegin;
112: if (st->transform) {
113: if (st->matmode == ST_MATMODE_COPY && nmat>2) {
114: nc = (nmat*(nmat+1))/2;
115: PetscCall(PetscMalloc1(nc,&coeffs));
116: /* Compute coeffs */
117: PetscCall(STCoeffs_Monomial(st,coeffs));
118: }
119: for (k=0;k<nmat-1;k++) PetscCall(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,PETSC_FALSE,&st->T[k]));
120: if (st->matmode == ST_MATMODE_COPY && nmat>2) PetscCall(PetscFree(coeffs));
121: if (st->nmat<=2) st->M = st->T[0];
122: }
123: PetscFunctionReturn(PETSC_SUCCESS);
124: }
126: SLEPC_EXTERN PetscErrorCode STCreate_Shift(ST st)
127: {
128: PetscFunctionBegin;
129: st->usesksp = PETSC_TRUE;
131: st->ops->apply = STApply_Generic;
132: st->ops->applytrans = STApplyTranspose_Generic;
133: st->ops->applyhermtrans = STApplyHermitianTranspose_Generic;
134: st->ops->backtransform = STBackTransform_Shift;
135: st->ops->setshift = STSetShift_Shift;
136: st->ops->getbilinearform = STGetBilinearForm_Default;
137: st->ops->setup = STSetUp_Shift;
138: st->ops->computeoperator = STComputeOperator_Shift;
139: st->ops->postsolve = STPostSolve_Shift;
140: st->ops->setdefaultksp = STSetDefaultKSP_Default;
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }