Actual source code: shift.c
slepc-3.22.2 2024-12-02
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: /*
18: Special STApply() for the BSE structured matrix
20: H = [ R C; -C^H -R^T ].
22: Assumes that H is a MATNEST and x,y are VECNEST.
24: Only the upper part of the product y=H*x is computed.
26: y1 = R*x1+C*x2
28: The bottom part of the input vector x2 should normally contain
29: either conj(x1) or -conj(x1).
30: The bottom part of the output vector y2 is not referenced.
31: */
32: static PetscErrorCode STApply_Shift_BSE(ST st,Vec x,Vec y)
33: {
34: Mat H,R,C;
35: Vec x1,x2,y1;
37: PetscFunctionBegin;
38: H = st->T[0];
39: PetscCall(MatNestGetSubMat(H,0,0,&R));
40: PetscCall(MatNestGetSubMat(H,0,1,&C));
41: PetscCall(VecNestGetSubVec(x,0,&x1));
42: PetscCall(VecNestGetSubVec(x,1,&x2));
43: PetscCall(VecNestGetSubVec(y,0,&y1));
44: PetscCall(MatMult(C,x2,y1));
45: PetscCall(MatMultAdd(R,x1,y1,y1));
46: PetscFunctionReturn(PETSC_SUCCESS);
47: }
49: static PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
50: {
51: PetscInt j;
53: PetscFunctionBegin;
54: for (j=0;j<n;j++) {
55: eigr[j] += st->sigma;
56: }
57: PetscFunctionReturn(PETSC_SUCCESS);
58: }
60: static PetscErrorCode STPostSolve_Shift(ST st)
61: {
62: PetscFunctionBegin;
63: if (st->matmode == ST_MATMODE_INPLACE) {
64: if (st->nmat>1) PetscCall(MatAXPY(st->A[0],st->sigma,st->A[1],st->str));
65: else PetscCall(MatShift(st->A[0],st->sigma));
66: st->Astate[0] = ((PetscObject)st->A[0])->state;
67: st->state = ST_STATE_INITIAL;
68: st->opready = PETSC_FALSE;
69: }
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
73: /*
74: Operator (shift):
75: Op P M
76: if nmat=1: A-sI NULL A-sI
77: if nmat=2: B^-1 (A-sB) B A-sB
78: */
79: static PetscErrorCode STComputeOperator_Shift(ST st)
80: {
81: PetscFunctionBegin;
82: st->usesksp = (st->nmat>1)? PETSC_TRUE: PETSC_FALSE;
83: PetscCall(PetscObjectReference((PetscObject)st->A[1]));
84: PetscCall(MatDestroy(&st->T[1]));
85: st->T[1] = st->A[1];
86: PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[0]));
87: if (st->nmat>1) PetscCall(PetscObjectReference((PetscObject)st->T[1]));
88: PetscCall(MatDestroy(&st->P));
89: st->P = (st->nmat>1)? st->T[1]: NULL;
90: st->M = st->T[0];
91: if (st->nmat>1 && st->Psplit) { /* build custom preconditioner from the split matrices */
92: PetscCall(MatDestroy(&st->Pmat));
93: PetscCall(PetscObjectReference((PetscObject)st->Psplit[1]));
94: st->Pmat = st->Psplit[1];
95: }
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: static PetscErrorCode STSetUp_Shift(ST st)
100: {
101: PetscInt k,nc,nmat=st->nmat;
102: PetscScalar *coeffs=NULL;
104: PetscFunctionBegin;
105: if (nmat>1) PetscCall(STSetWorkVecs(st,1));
106: if (nmat>2) { /* set-up matrices for polynomial eigenproblems */
107: if (st->transform) {
108: nc = (nmat*(nmat+1))/2;
109: PetscCall(PetscMalloc1(nc,&coeffs));
110: /* Compute coeffs */
111: PetscCall(STCoeffs_Monomial(st,coeffs));
112: /* T[n] = A_n */
113: k = nmat-1;
114: PetscCall(PetscObjectReference((PetscObject)st->A[k]));
115: PetscCall(MatDestroy(&st->T[k]));
116: st->T[k] = st->A[k];
117: 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]));
118: PetscCall(PetscFree(coeffs));
119: PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
120: PetscCall(MatDestroy(&st->P));
121: st->P = st->T[nmat-1];
122: if (st->Psplit) { /* build custom preconditioner from the split matrices */
123: PetscCall(STMatMAXPY_Private(st,st->sigma,0.0,nmat-1,coeffs?coeffs:NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
124: }
125: PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
126: } else {
127: for (k=0;k<nmat;k++) {
128: PetscCall(PetscObjectReference((PetscObject)st->A[k]));
129: PetscCall(MatDestroy(&st->T[k]));
130: st->T[k] = st->A[k];
131: }
132: }
133: }
134: if (st->P) PetscCall(KSPSetUp(st->ksp));
135: if (st->structured) st->ops->apply = STApply_Shift_BSE;
136: PetscFunctionReturn(PETSC_SUCCESS);
137: }
139: static PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
140: {
141: PetscInt k,nc,nmat=PetscMax(st->nmat,2);
142: PetscScalar *coeffs=NULL;
144: PetscFunctionBegin;
145: if (st->transform) {
146: if (st->matmode == ST_MATMODE_COPY && nmat>2) {
147: nc = (nmat*(nmat+1))/2;
148: PetscCall(PetscMalloc1(nc,&coeffs));
149: /* Compute coeffs */
150: PetscCall(STCoeffs_Monomial(st,coeffs));
151: }
152: 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]));
153: if (st->matmode == ST_MATMODE_COPY && nmat>2) PetscCall(PetscFree(coeffs));
154: if (st->nmat<=2) st->M = st->T[0];
155: }
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: SLEPC_EXTERN PetscErrorCode STCreate_Shift(ST st)
160: {
161: PetscFunctionBegin;
162: st->usesksp = PETSC_TRUE;
164: st->ops->apply = STApply_Generic;
165: st->ops->applytrans = STApplyTranspose_Generic;
166: st->ops->applyhermtrans = STApplyHermitianTranspose_Generic;
167: st->ops->backtransform = STBackTransform_Shift;
168: st->ops->setshift = STSetShift_Shift;
169: st->ops->getbilinearform = STGetBilinearForm_Default;
170: st->ops->setup = STSetUp_Shift;
171: st->ops->computeoperator = STComputeOperator_Shift;
172: st->ops->postsolve = STPostSolve_Shift;
173: st->ops->setdefaultksp = STSetDefaultKSP_Default;
174: PetscFunctionReturn(PETSC_SUCCESS);
175: }